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

AlbPUMI_FMDBMeshStruct.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 
00007 #include <iostream>
00008 
00009 #include <boost/mpi/collectives/all_gather.hpp>
00010 #include <boost/mpi/collectives/all_reduce.hpp>
00011 #include "AlbPUMI_FMDBMeshStruct.hpp"
00012 
00013 #include "Teuchos_VerboseObject.hpp"
00014 #include "Albany_Utils.hpp"
00015 
00016 #include "Teuchos_TwoDArray.hpp"
00017 #include <Shards_BasicTopologies.hpp>
00018 
00019 #include <apfSTK.h>
00020 #include <apfShape.h>
00021 #include <ma.h>
00022 
00023 class SizeFunction : public ma::IsotropicFunction {
00024   public:
00025     SizeFunction(double s) {size = s;}
00026     double getValue(ma::Entity*) {return size;}
00027   private:
00028     double size;
00029 };
00030 
00031 AlbPUMI::FMDBMeshStruct::FMDBMeshStruct(
00032           const Teuchos::RCP<Teuchos::ParameterList>& params,
00033       const Teuchos::RCP<const Epetra_Comm>& comm) :
00034   out(Teuchos::VerboseObjectBase::getDefaultOStream()),
00035   apfMesh(0)
00036 {
00037   // fmdb skips mpi initialization if it's already initialized
00038   SCUTIL_Init(Albany::getMpiCommFromEpetraComm(*comm));
00039   params->validateParameters(*getValidDiscretizationParameters(),0);
00040 
00041   std::string mesh_file = params->get<std::string>("FMDB Input File Name");
00042   outputFileName = params->get<std::string>("FMDB Output File Name", "");
00043   outputInterval = params->get<int>("FMDB Write Interval", 1); // write every time step default
00044   if (params->get<bool>("Call serial global partition"))
00045     useDistributedMesh=false;
00046   else
00047     useDistributedMesh=true;
00048 
00049   compositeTet = params->get<bool>("Use Composite Tet 10", false);
00050 
00051   // create a model and load
00052   model = NULL; // default is no model
00053 
00054   PUMI_Geom_RegisterMesh();
00055 
00056   if(params->isParameter("Acis Model Input File Name")){ // User has an Acis model
00057 
00058     std::string model_file = params->get<std::string>("Acis Model Input File Name");
00059     PUMI_Geom_RegisterAcis();
00060     PUMI_Geom_LoadFromFile(model, model_file.c_str());
00061   }
00062 
00063   if(params->isParameter("Parasolid Model Input File Name")){ // User has a Parasolid model
00064 
00065     std::string model_file = params->get<std::string>("Parasolid Model Input File Name");
00066     PUMI_Geom_RegisterParasolid();
00067     PUMI_Geom_LoadFromFile(model, model_file.c_str());
00068   }
00069 
00070   if(params->isParameter("Mesh Model Input File Name")){ // User has a meshModel model
00071 
00072     std::string model_file = params->get<std::string>("Mesh Model Input File Name");
00073     PUMI_Geom_LoadFromFile(model, model_file.c_str());
00074   }
00075 
00076   TEUCHOS_TEST_FOR_EXCEPTION(model==NULL,std::logic_error,"FMDBMeshStruct: no model" << std::endl);
00077 
00078   FMDB_Mesh_Create (model, mesh);
00079 
00080   int rc = FMDB_Mesh_LoadFromFile (mesh, &mesh_file[0], useDistributedMesh);
00081   TEUCHOS_TEST_FOR_EXCEPTION(rc,std::logic_error,
00082       "FAILED MESH LOADING - check mesh file or number of input files" << std::endl)
00083 
00084   FMDB_Mesh_DspSize(mesh);
00085 
00086   if(params->isParameter("Element Block Associations")){ // User has specified associations in the input file
00087 
00088     // Get element block associations from input file
00089     Teuchos::TwoDArray< std::string > EBAssociations;
00090 
00091     EBAssociations = params->get<Teuchos::TwoDArray<std::string> >("Element Block Associations");
00092 
00093     TEUCHOS_TEST_FOR_EXCEPTION( !(2 == EBAssociations.getNumRows()),
00094             Teuchos::Exceptions::InvalidParameter,
00095             "Error in specifying element block associations in input file" );
00096 
00097     int nEBAssoc = EBAssociations.getNumCols();
00098 
00099     for(size_t eb = 0; eb < nEBAssoc; eb++){
00100       *out << "Element block \"" <<  EBAssociations(1, eb).c_str() << "\" matches mesh region : "
00101            << EBAssociations(0, eb).c_str() << std::endl;
00102     }
00103 
00104     pumi::GRIter gr_iter = pumi::GM_regionIter(model);
00105     pGeomEnt geom_rgn;
00106     while ((geom_rgn = pumi::GRIter_next(gr_iter)) != NULL)
00107     {
00108       for(size_t eblock = 0; eblock < nEBAssoc; eblock++){
00109         if (GEN_tag(geom_rgn) == atoi(EBAssociations(0, eblock).c_str()))
00110           PUMI_Exodus_CreateElemBlk(geom_rgn, EBAssociations(1, eblock).c_str());
00111       }
00112     }
00113     pumi::GRIter_delete(gr_iter);
00114   }
00115 
00116 
00117   if(params->isParameter("Node Set Associations")){ // User has specified associations in the input file
00118 
00119     // Get node set associations from input file
00120     Teuchos::TwoDArray< std::string > NSAssociations;
00121 
00122     NSAssociations = params->get<Teuchos::TwoDArray<std::string> >("Node Set Associations");
00123 
00124     TEUCHOS_TEST_FOR_EXCEPTION( !(2 == NSAssociations.getNumRows()),
00125         Teuchos::Exceptions::InvalidParameter,
00126         "Error in specifying node set associations in input file" );
00127 
00128     int nNSAssoc = NSAssociations.getNumCols();
00129 
00130     for(size_t ns = 0; ns < nNSAssoc; ns++){
00131       *out << "Node set \"" << NSAssociations(1, ns).c_str() << "\" matches geometric face : "
00132         << NSAssociations(0, ns).c_str() << std::endl;
00133     }
00134 
00135     pumi::GFIter gf_iter=pumi::GM_faceIter(model);
00136     pGeomEnt geom_face;
00137     while ((geom_face=pumi::GFIter_next(gf_iter)) != NULL)
00138     {
00139       for(size_t ns = 0; ns < nNSAssoc; ns++){
00140         if (GEN_tag(geom_face) == atoi(NSAssociations(0, ns).c_str())){
00141           PUMI_Exodus_CreateNodeSet(geom_face, NSAssociations(1, ns).c_str());
00142         }
00143       }
00144     }
00145     pumi::GFIter_delete(gf_iter);
00146 
00147   }
00148 
00149   if(params->isParameter("Edge Node Set Associations")){ // User has specified associations in the input file
00150 
00151     // Get node set associations from input file
00152     Teuchos::TwoDArray< std::string > EdgeNSAssociations;
00153 
00154     EdgeNSAssociations = params->get<Teuchos::TwoDArray<std::string> >("Edge Node Set Associations");
00155 
00156     TEUCHOS_TEST_FOR_EXCEPTION( !(2 == EdgeNSAssociations.getNumRows()),
00157         Teuchos::Exceptions::InvalidParameter,
00158         "Error in specifying node set associations in input file" );
00159 
00160     int nEdgeNSAssoc = EdgeNSAssociations.getNumCols();
00161 
00162     for(size_t ns = 0; ns < nEdgeNSAssoc; ns++){
00163       *out << "Node set \"" << EdgeNSAssociations(1, ns).c_str() << "\" matches geometric edge : "
00164         << EdgeNSAssociations(0, ns).c_str() << std::endl;
00165     }
00166 
00167     pumi::GEIter ge_iter=pumi::GM_edgeIter(model);
00168     pGeomEnt geom_edge;
00169     while ((geom_edge=pumi::GEIter_next(ge_iter)) != NULL)
00170     {
00171       for(size_t ns = 0; ns < nEdgeNSAssoc; ns++){
00172         if (GEN_tag(geom_edge) == atoi(EdgeNSAssociations(0, ns).c_str())){
00173           PUMI_Exodus_CreateNodeSet(geom_edge, EdgeNSAssociations(1, ns).c_str());
00174         }
00175       }
00176     }
00177     pumi::GEIter_delete(ge_iter);
00178   }
00179 
00180   if(params->isParameter("Vertex Node Set Associations")){ // User has specified associations in the input file
00181 
00182     // Get node set associations from input file
00183     Teuchos::TwoDArray< std::string > VertexNSAssociations;
00184 
00185     VertexNSAssociations = params->get<Teuchos::TwoDArray<std::string> >("Vertex Node Set Associations");
00186 
00187     TEUCHOS_TEST_FOR_EXCEPTION( !(2 == VertexNSAssociations.getNumRows()),
00188         Teuchos::Exceptions::InvalidParameter,
00189         "Error in specifying node set associations in input file" );
00190 
00191     int nVertexNSAssoc = VertexNSAssociations.getNumCols();
00192 
00193     for(size_t ns = 0; ns < nVertexNSAssoc; ns++){
00194       *out << "Node set \"" << VertexNSAssociations(1, ns).c_str() << "\" matches geometric vertex : "
00195         << VertexNSAssociations(0, ns).c_str() << std::endl;
00196     }
00197 
00198     pumi::GVIter gv_iter=pumi::GM_vertexIter(model);
00199     pGeomEnt geom_vertex;
00200     while ((geom_vertex=pumi::GVIter_next(gv_iter)) != NULL)
00201     {
00202       for(size_t ns = 0; ns < nVertexNSAssoc; ns++){
00203         if (GEN_tag(geom_vertex) == atoi(VertexNSAssociations(0, ns).c_str())){
00204           PUMI_Exodus_CreateNodeSet(geom_vertex, VertexNSAssociations(1, ns).c_str());
00205         }
00206       }
00207     }
00208     pumi::GVIter_delete(gv_iter);
00209   }
00210 
00211   if(params->isParameter("Side Set Associations")){ // User has specified associations in the input file
00212 
00213     // Get side set block associations from input file
00214     Teuchos::TwoDArray< std::string > SSAssociations;
00215 
00216     SSAssociations = params->get<Teuchos::TwoDArray<std::string> >("Side Set Associations");
00217 
00218     TEUCHOS_TEST_FOR_EXCEPTION( !(2 == SSAssociations.getNumRows()),
00219             Teuchos::Exceptions::InvalidParameter,
00220             "Error in specifying side set associations in input file" );
00221 
00222     int nSSAssoc = SSAssociations.getNumCols();
00223 
00224     for(size_t ss = 0; ss < nSSAssoc; ss++){
00225       *out << "Side set \"" << SSAssociations(1, ss).c_str() << "\" matches geometric face : "
00226            << SSAssociations(0, ss).c_str() << std::endl;
00227     }
00228 
00229 
00230     pumi::GFIter gf_iter=pumi::GM_faceIter(model);
00231     pGeomEnt geom_face;
00232     while ((geom_face = pumi::GFIter_next(gf_iter)) != NULL)
00233     {
00234       for(size_t ss = 0; ss < nSSAssoc; ss++){
00235         if (GEN_tag(geom_face) == atoi(SSAssociations(0, ss).c_str()))
00236           PUMI_Exodus_CreateSideSet(geom_face, SSAssociations(1, ss).c_str());
00237       }
00238     }
00239     pumi::GFIter_delete(gf_iter);
00240   }
00241 
00242   apfMesh = apf::createMesh(mesh);
00243   bool isQuadMesh = params->get<bool>("2nd Order Mesh",false);
00244   if (isQuadMesh)
00245     changeMeshShape(apfMesh,apf::getLagrange(2));
00246 
00247   // Resize mesh after input if indicated in the input file
00248   if(params->isParameter("Resize Input Mesh Element Size")){ // User has indicated a desired element size in input file
00249       SizeFunction sizeFunction(params->get<double>("Resize Input Mesh Element Size", 0.1));
00250       int num_iters = params->get<int>("Max Number of Mesh Adapt Iterations", 1);
00251       ma::Input* input = ma::configure(apfMesh,&sizeFunction);
00252       input->maximumIterations = num_iters;
00253       input->shouldSnap = false;
00254       ma::adapt(input);
00255       FMDB_Mesh_DspSize(mesh);
00256   }
00257 
00258   //get mesh dim
00259   FMDB_Mesh_GetDim(mesh, &numDim);
00260 
00261   std::vector<pElemBlk> elem_blocks;
00262   PUMI_Exodus_GetElemBlk(mesh, elem_blocks);
00263 
00264   // Build a map to get the EB name given the index
00265 
00266   int numEB = elem_blocks.size(), EB_size;
00267   *out <<"["<<SCUTIL_CommRank()<< "] Found : " << numEB << " element blocks." << std::endl;
00268   std::vector<int> el_blocks;
00269 
00270   for (int eb=0; eb < numEB; eb++){
00271     std::string EB_name;
00272     PUMI_ElemBlk_GetName(elem_blocks[eb], EB_name);
00273     this->ebNameToIndex[EB_name] = eb;
00274     PUMI_ElemBlk_GetSize(mesh, elem_blocks[eb], &EB_size);
00275     el_blocks.push_back(EB_size);
00276   }
00277 
00278   // Set defaults for cubature and workset size, overridden in input file
00279 
00280   cubatureDegree = params->get("Cubature Degree", 1);
00281   int worksetSizeMax = params->get("Workset Size", 50);
00282   interleavedOrdering = params->get("Interleaved Ordering",true);
00283   allElementBlocksHaveSamePhysics = true;
00284   hasRestartSolution = false;
00285 
00286   // No history available by default
00287   solutionFieldHistoryDepth = 0;
00288 
00289   // This is typical, can be resized for multiple material problems
00290   meshSpecs.resize(1);
00291 
00292   // Get number of elements per element block
00293   // in calculating an upper bound on the worksetSize.
00294 
00295   int ebSizeMax =  *std::max_element(el_blocks.begin(), el_blocks.end());
00296   worksetSize = computeWorksetSize(worksetSizeMax, ebSizeMax);
00297 
00298   *out <<"["<<SCUTIL_CommRank()<< "] Workset size is: " << worksetSize << std::endl;
00299 
00300   // Node sets
00301   std::vector<pNodeSet> node_sets;
00302   PUMI_Exodus_GetNodeSet(mesh, node_sets);
00303 
00304   for(int ns = 0; ns < node_sets.size(); ns++)
00305   {
00306     std::string NS_name;
00307     PUMI_NodeSet_GetName(node_sets[ns], NS_name);
00308     nsNames.push_back(NS_name);
00309   }
00310 
00311   // Side sets
00312   std::vector<pSideSet> side_sets;
00313   PUMI_Exodus_GetSideSet(mesh, side_sets);
00314 
00315   int status;
00316 
00317   for(int ss = 0; ss < side_sets.size(); ss++)
00318   {
00319     std::string SS_name;
00320     PUMI_SideSet_GetName(side_sets[ss], SS_name);
00321     ssNames.push_back(SS_name);
00322   }
00323 
00324   // Construct MeshSpecsStruct
00325   std::vector<pMeshEnt> elements;
00326   const CellTopologyData* ctd = getCellTopology(apfMesh);
00327   if (!params->get("Separate Evaluators by Element Block",false))
00328   {
00329     // get elements in the first element block
00330     std::string EB_name;
00331     PUMI_ElemBlk_GetName(elem_blocks[0], EB_name);
00332     this->meshSpecs[0] = Teuchos::rcp(
00333         new Albany::MeshSpecsStruct(
00334           *ctd, numDim, cubatureDegree,
00335           nsNames, ssNames, worksetSize, EB_name,
00336           this->ebNameToIndex, this->interleavedOrdering));
00337   }
00338   else
00339   {
00340     *out <<"["<<SCUTIL_CommRank()<< "] MULTIPLE Elem Block in FMDB: DO worksetSize[eb] max?? " << std::endl;
00341     this->allElementBlocksHaveSamePhysics=false;
00342     this->meshSpecs.resize(numEB);
00343     int eb_size;
00344     std::string eb_name;
00345     for (int eb=0; eb<numEB; eb++)
00346     {
00347       std::string EB_name;
00348       PUMI_ElemBlk_GetName(elem_blocks[eb], EB_name);
00349       this->meshSpecs[eb] = Teuchos::rcp(new Albany::MeshSpecsStruct(*ctd, numDim, cubatureDegree,
00350                                               nsNames, ssNames, worksetSize, EB_name,
00351                                               this->ebNameToIndex, this->interleavedOrdering));
00352     } // for
00353   } // else
00354 
00355 }
00356 
00357 AlbPUMI::FMDBMeshStruct::~FMDBMeshStruct()
00358 {
00359   apfMesh->destroyNative();
00360   apf::destroyMesh(apfMesh);
00361   ParUtil::Instance()->Finalize(0); // skip MPI_finalize
00362 }
00363 
00364 void
00365 AlbPUMI::FMDBMeshStruct::setFieldAndBulkData(
00366                   const Teuchos::RCP<const Epetra_Comm>& comm,
00367                   const Teuchos::RCP<Teuchos::ParameterList>& params,
00368                   const unsigned int neq_,
00369                   const Albany::AbstractFieldContainer::FieldContainerRequirements& req,
00370                   const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00371                   const unsigned int worksetSize_)
00372 {
00373 
00374   using Albany::StateStruct;
00375 
00376   // Set the number of equation present per node. Needed by AlbPUMI_FMDBDiscretization.
00377 
00378   neq = neq_;
00379 
00380   Teuchos::Array<std::string> defaultLayout;
00381   solVectorLayout = 
00382     params->get<Teuchos::Array<std::string> >("Solution Vector Components", defaultLayout);
00383 
00384   if (solVectorLayout.size() == 0) { 
00385     int valueType;
00386     if (neq==1)
00387       valueType = apf::SCALAR;
00388     else if (neq==3)
00389       valueType = apf::VECTOR;
00390     else
00391     {
00392       assert(neq==9);
00393       valueType = apf::MATRIX;
00394     }
00395     apf::createFieldOn(apfMesh,"residual",valueType);
00396     apf::createFieldOn(apfMesh,"solution",valueType);
00397   }
00398   else 
00399     splitFields(solVectorLayout);
00400 
00401   solutionInitialized = false;
00402   residualInitialized = false;
00403 
00404   // Code to parse the vector of StateStructs and save the information
00405 
00406   // dim[0] is the number of cells in this workset
00407   // dim[1] is the number of QP per cell
00408   // dim[2] is the number of dimensions of the field
00409   // dim[3] is the number of dimensions of the field
00410 
00411   std::set<std::string> nameSet;
00412 
00413   for (std::size_t i=0; i<sis->size(); i++) {
00414     StateStruct& st = *((*sis)[i]);
00415 
00416     if ( ! nameSet.insert(st.name).second)
00417       continue; //ignore duplicates
00418 
00419     std::vector<int>& dim = st.dim;
00420 
00421     if(st.entity == StateStruct::NodalData) { // Data at the node points
00422 
00423        const Teuchos::RCP<Albany::NodeFieldContainer>& nodeContainer
00424                = sis->getNodalDataBlock()->getNodeContainer();
00425 
00426         (*nodeContainer)[st.name] = AlbPUMI::buildPUMINodeField(st.name, dim, st.output);
00427 
00428     }
00429 
00430     // qpscalars
00431 
00432     else if (dim.size() == 2){
00433       if(st.entity == StateStruct::QuadPoint || st.entity == StateStruct::ElemNode) {
00434 
00435         qpscalar_states.push_back(Teuchos::rcp(new QPData<double, 2>(st.name, dim, st.output)));
00436 
00437         if ( ! PCU_Comm_Self())
00438           std::cout << "AlbPUMI::FMDBMeshStruct qps field name " << st.name
00439             << " size : " << dim[1] << std::endl;
00440       }
00441     }
00442 
00443     // qpvectors
00444 
00445     else if (dim.size() == 3){
00446       if(st.entity == StateStruct::QuadPoint || st.entity == StateStruct::ElemNode) {
00447 
00448         qpvector_states.push_back(Teuchos::rcp(new QPData<double, 3>(st.name, dim, st.output)));
00449 
00450         if ( ! PCU_Comm_Self())
00451           std::cout << "AlbPUMI::FMDBMeshStruct qpv field name " << st.name
00452             << " dim[1] : " << dim[1] << " dim[2] : " << dim[2] << std::endl;
00453       }
00454     }
00455 
00456     // qptensors
00457 
00458     else if (dim.size() == 4){
00459       if(st.entity == StateStruct::QuadPoint || st.entity == StateStruct::ElemNode) {
00460 
00461         qptensor_states.push_back(Teuchos::rcp(new QPData<double, 4>(st.name, dim, st.output)));
00462 
00463         if ( ! PCU_Comm_Self())
00464           std::cout << "AlbPUMI::FMDBMeshStruct qpt field name " << st.name
00465             << " dim[1] : " << dim[1] << " dim[2] : " << dim[2] << " dim[3] : " << dim[3] << std::endl;
00466       }
00467     }
00468 
00469     // just a scalar number
00470 
00471     else if ( dim.size() == 1 && st.entity == Albany::StateStruct::WorksetValue) {
00472       // dim not used or accessed here
00473       scalarValue_states.push_back(Teuchos::rcp(new QPData<double, 1>(st.name, dim, st.output)));
00474     }
00475 
00476     // anything else is an error!
00477 
00478     else TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "dim.size() < 2 || dim.size()>4 || " <<
00479          "st.entity != Albany::StateStruct::QuadPoint || " <<
00480          "st.entity != Albany::StateStruct::ElemNode || " <<
00481          "st.entity != Albany::StateStruct::NodalData" << std::endl);
00482 
00483   }
00484 
00485 }
00486 
00487 void
00488 AlbPUMI::FMDBMeshStruct::splitFields(Teuchos::Array<std::string> fieldLayout)
00489 { // user is breaking up or renaming solution & residual fields
00490 
00491   TEUCHOS_TEST_FOR_EXCEPTION((fieldLayout.size() % 2), std::logic_error,
00492       "Error in input file: specification of solution vector layout is incorrect\n");
00493 
00494   int valueType;
00495 
00496   for (std::size_t i=0; i < fieldLayout.size(); i+=2) {
00497 
00498     if (fieldLayout[i+1] == "S")
00499       valueType = apf::SCALAR;
00500     else if (fieldLayout[i+1] == "V")
00501       valueType = apf::VECTOR;
00502     else
00503       TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
00504           "Error in input file: specification of solution vector layout is incorrect\n");
00505 
00506     apf::createFieldOn(apfMesh,fieldLayout[i].c_str(),valueType);
00507     apf::createFieldOn(apfMesh,fieldLayout[i].append("Res").c_str(),valueType);
00508   }
00509 
00510 }
00511 
00512 Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >&
00513 AlbPUMI::FMDBMeshStruct::getMeshSpecs()
00514 {
00515   TEUCHOS_TEST_FOR_EXCEPTION(meshSpecs==Teuchos::null,
00516        std::logic_error,
00517        "meshSpecs accessed, but it has not been constructed" << std::endl);
00518   return meshSpecs;
00519 }
00520 
00521 Albany::AbstractMeshStruct::msType
00522 AlbPUMI::FMDBMeshStruct::meshSpecsType()
00523 {
00524 
00525   std::string str = outputFileName;
00526   size_t found = str.find("vtk");
00527 
00528   if(found != std::string::npos){
00529 
00530     return FMDB_VTK_MS;
00531 
00532   }
00533 
00534   found = str.find("exo");
00535   if(found != std::string::npos){
00536 
00537     return FMDB_EXODUS_MS;
00538 
00539   }
00540 
00541   TEUCHOS_TEST_FOR_EXCEPTION(true,
00542        std::logic_error,
00543        "Unrecognized output file extension given in the input file" << std::endl);
00544 
00545 }
00546 
00547 int AlbPUMI::FMDBMeshStruct::computeWorksetSize(const int worksetSizeMax,
00548                                                      const int ebSizeMax) const
00549 {
00550   // Resize workset size down to maximum number in an element block
00551   if (worksetSizeMax > ebSizeMax || worksetSizeMax < 1) return ebSizeMax;
00552   else {
00553      // compute numWorksets, and shrink workset size to minimize padding
00554      const int numWorksets = 1 + (ebSizeMax-1) / worksetSizeMax;
00555      return (1 + (ebSizeMax-1) /  numWorksets);
00556   }
00557 }
00558 
00559 void
00560 AlbPUMI::FMDBMeshStruct::loadSolutionFieldHistory(int step)
00561 {
00562   TEUCHOS_TEST_FOR_EXCEPT(step < 0 || step >= solutionFieldHistoryDepth);
00563 }
00564 
00565 void AlbPUMI::FMDBMeshStruct::setupMeshBlkInfo()
00566 {
00567 
00568    int nBlocks = this->meshSpecs.size();
00569 
00570    for(int i = 0; i < nBlocks; i++){
00571 
00572       const Albany::MeshSpecsStruct &ms = *meshSpecs[i];
00573 
00574       meshDynamicData[i] = Teuchos::rcp(new Albany::CellSpecs(ms.ctd, ms.worksetSize, ms.cubatureDegree,
00575                       numDim, neq, 0, useCompositeTet()));
00576 
00577    }
00578 
00579 }
00580 
00581 Teuchos::RCP<const Teuchos::ParameterList>
00582 AlbPUMI::FMDBMeshStruct::getValidDiscretizationParameters() const
00583 {
00584 
00585   Teuchos::RCP<Teuchos::ParameterList> validPL
00586      = rcp(new Teuchos::ParameterList("Valid FMDBParams"));
00587 
00588   validPL->set<std::string>("FMDB Solution Name", "",
00589       "Name of solution output vector written to output file");
00590   validPL->set<std::string>("FMDB Residual Name", "",
00591       "Name of residual output vector written to output file");
00592   validPL->set<int>("FMDB Write Interval", 3, "Step interval to write solution data to output file");
00593   validPL->set<std::string>("Method", "",
00594     "The discretization method, parsed in the Discretization Factory");
00595   validPL->set<int>("Cubature Degree", 1, "Integration order sent to Intrepid");
00596   validPL->set<int>("Workset Size", 50, "Upper bound on workset (bucket) size");
00597   validPL->set<bool>("Interleaved Ordering", true, "Flag for interleaved or blocked unknown ordering");
00598   validPL->set<bool>("Separate Evaluators by Element Block", false,
00599                      "Flag for different evaluation trees for each Element Block");
00600   Teuchos::Array<std::string> defaultFields;
00601   validPL->set<Teuchos::Array<std::string> >("Restart Fields", defaultFields,
00602       "Fields to pick up from the restart file when restarting");
00603   validPL->set<Teuchos::Array<std::string> >("Solution Vector Components", defaultFields,
00604       "Names and layouts of solution vector components");
00605   validPL->set<bool>("2nd Order Mesh", false, "Flag to indicate 2nd order Lagrange shape functions");
00606   
00607   validPL->set<std::string>("FMDB Input File Name", "", "File Name For FMDB Mesh Input");
00608   validPL->set<std::string>("FMDB Output File Name", "", "File Name For FMDB Mesh Output");
00609 
00610   validPL->set<std::string>("Acis Model Input File Name", "", "File Name For ACIS Model Input");
00611   validPL->set<std::string>("Parasolid Model Input File Name", "", "File Name For PARASOLID Model Input");
00612   validPL->set<std::string>("Mesh Model Input File Name", "", "File Name for meshModel Input");
00613 
00614   validPL->set<bool>("Periodic BC", false, "Flag to indicate periodic a mesh");
00615   validPL->set<int>("Restart Index", 1, "Exodus time index to read for initial guess/condition.");
00616   validPL->set<double>("Restart Time", 1.0, "Exodus solution time to read for initial guess/condition.");
00617   validPL->set<bool>("Use Serial Mesh", false, "Read in a single mesh on PE 0 and rebalance");
00618 
00619   validPL->set<int>("Number of parts", 1, "Number of parts");
00620   validPL->set<int>("Number of migrations", 0, "Number of migrations");
00621   validPL->set<int>("Number of individual migrations", 0, "Number of individual migrations");
00622   validPL->set<double>("Imbalance tolerance", 1.03, "Imbalance tolerance");
00623   validPL->set<bool>("Construct pset", false, "Construct pset");
00624   validPL->set<bool>("Call serial global partition", false, "Call serial global partition");
00625 
00626   // Parameters to refine the mesh after input
00627   validPL->set<double>("Resize Input Mesh Element Size", 1.0, "Resize mesh element to this size at input");
00628   validPL->set<int>("Max Number of Resize Iterations", 0, "Max number of iteration sweeps to use during initial element resize");
00629 
00630   validPL->set<std::string>("LB Method", "", "Method used to load balance mesh (default \"ParMETIS\")");
00631   validPL->set<std::string>("LB Approach", "", "Approach used to load balance mesh (default \"PartKway\")");
00632 
00633   Teuchos::TwoDArray<std::string> defaultData;
00634   validPL->set<Teuchos::TwoDArray<std::string> >("Element Block Associations", defaultData,
00635       "Association between region ID and element block string");
00636   validPL->set<Teuchos::TwoDArray<std::string> >("Node Set Associations", defaultData,
00637       "Association between geometric face ID and node set string");
00638   validPL->set<Teuchos::TwoDArray<std::string> >("Edge Node Set Associations", defaultData,
00639       "Association between geometric edge ID and node set string");
00640   validPL->set<Teuchos::TwoDArray<std::string> >("Vertex Node Set Associations", defaultData,
00641       "Association between geometric edge ID and node set string");
00642   validPL->set<Teuchos::TwoDArray<std::string> >("Side Set Associations", defaultData,
00643       "Association between face ID and side set string");
00644 
00645   return validPL;
00646 }
00647 

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