FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

fieldview9Reader.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Description
00025     Reader module for Fieldview9 to read OpenFOAM mesh and data.
00026 
00027     Creates new 'fvbin' type executable which needs to be installed in place
00028     of bin/fvbin.
00029 
00030     Implements a reader for combined mesh&results on an unstructured mesh.
00031 
00032     See Fieldview Release 9 Reference Manual and coding in user/ directory
00033     of the Fieldview release.
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #include <finiteVolume/fvCFD.H>
00038 #include <OpenFOAM/IOobjectList.H>
00039 #include <OpenFOAM/GeometricField.H>
00040 #include <OpenFOAM/pointFields.H>
00041 #include <finiteVolume/volPointInterpolation.H>
00042 #include "readerDatabase.H"
00043 #include <OpenFOAM/wallPolyPatch.H>
00044 #include <OpenFOAM/ListOps.H>
00045 #include <OpenFOAM/cellModeller.H>
00046 
00047 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00048 
00049 extern "C"
00050 {
00051 // Define various Fieldview constants and prototypes
00052 
00053 #include "fv_reader_tags.h"
00054 
00055 static const int FVHEX     = 2;
00056 static const int FVPRISM   = 3;
00057 static const int FVPYRAMID = 4;
00058 static const int FVTET     = 1;
00059 static const int ITYPE     = 1;
00060 
00061 unsigned int fv_encode_elem_header(int elem_type, int wall_info[]);
00062 void time_step_get_value(float*, int, int*, float*, int*);
00063 void fv_error_msg(const char*, const char*);
00064 
00065 void reg_single_unstruct_reader
00066 (
00067     char *,
00068     void
00069     (
00070         char*, int*, int*, int*, int*, int*, int*,
00071         int[], int*, char[][80], int[], int*, char[][80], int*
00072     ),
00073     void
00074     (
00075         int*, int*, int*, float[], int*, float[], int*
00076     )
00077 );
00078 
00079 int create_tet(const int, const int[8], const int[]);
00080 int create_pyramid(const int, const int[8], const int[]);
00081 int create_prism(const int, const int[8], const int[]);
00082 int create_hex(const int, const int[8], const int[]);
00083 
00084 typedef unsigned char uChar;
00085 extern uChar create_bndry_face_buffered
00086 (
00087     int bndry_type,
00088     int num_verts,
00089     int verts[],
00090     int *normals_flags,
00091     int num_grid_nodes
00092 );
00093 
00094 /*
00095  * just define empty readers here for ease of linking.
00096  * Comment out if you have doubly defined linking error on this symbol
00097  */
00098 void ftn_register_data_readers()
00099 {}
00100 
00101 /*
00102  * just define empty readers here for ease of linking.
00103  * Comment out if you have doubly defined linking error on this symbol
00104  */
00105 void ftn_register_functions()
00106 {}
00107 
00108 /*
00109  * just define empty readers here for ease of linking.
00110  * Comment out if you have doubly defined linking error on this symbol
00111  */
00112 void register_functions()
00113 {}
00114 
00115 
00116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00117 
00118 
00119 //
00120 // Storage for all Foam state (mainly database & mesh)
00121 //
00122 static readerDatabase db_;
00123 
00124 
00125 // Write fv error message.
00126 static void errorMsg(const string& msg)
00127 {
00128     fv_error_msg("Foam Reader", msg.c_str());
00129 }
00130 
00131 
00132 // Simple check if directory is valid case directory.
00133 static bool validCase(const fileName& rootAndCase)
00134 {
00135     //if (isDir(rootAndCase/"system") && isDir(rootAndCase/"constant"))
00136     if (isDir(rootAndCase/"constant"))
00137     {
00138         return true;
00139     }
00140     else
00141     {
00142         return false;
00143     }
00144 }
00145 
00146 
00147 // Check whether case has topo changes by looking back from last time
00148 // to first directory with polyMesh/cells.
00149 static bool hasTopoChange(const instantList& times)
00150 {
00151     label lastIndex = times.size()-1;
00152 
00153     const Time& runTime = db_.runTime();
00154 
00155     // Only set time; do not update mesh.
00156     runTime.setTime(times[lastIndex], lastIndex);
00157 
00158     fileName facesInst(runTime.findInstance(polyMesh::meshSubDir, "faces"));
00159 
00160     // See if cellInst is constant directory. Note extra .name() is for
00161     // parallel cases when runTime.constant is '../constant'
00162     if (facesInst != runTime.constant().name())
00163     {
00164         Info<< "Found cells file in " << facesInst << " so case has "
00165             << "topological changes" << endl;
00166 
00167         return true;
00168     }
00169     else
00170     {
00171         Info<< "Found cells file in " << facesInst << " so case has "
00172             << "no topological changes" << endl;
00173 
00174         return false;
00175     }
00176 }
00177 
00178 
00179 static bool selectTime(const instantList& times, int* iret)
00180 {
00181     List<float> fvTimes(2*times.size());
00182 
00183     forAll(times, timeI)
00184     {
00185         fvTimes[2*timeI]   = float(timeI);
00186         fvTimes[2*timeI+1] = float(times[timeI].value());
00187     }
00188 
00189     int istep;
00190 
00191     float time;
00192 
00193     *iret=0;
00194 
00195     time_step_get_value(fvTimes.begin(), times.size(), &istep, &time, iret);
00196 
00197     if (*iret == -5)
00198     {
00199         errorMsg("Out of memory.");
00200 
00201         return false;
00202     }
00203     if (*iret == -15)
00204     {
00205         // Cancel action.
00206         return false;
00207     }
00208     if (*iret != 0)
00209     {
00210         errorMsg("Unspecified error.");
00211 
00212         return false;
00213     }
00214     Info<< "Selected timeStep:" << istep << "  time:" << scalar(time) << endl;
00215 
00216     // Set time and load mesh.
00217     db_.setTime(times[istep], istep);
00218 
00219     return true;
00220 }
00221 
00222 
00223 // Gets (names of) all fields in all timesteps.
00224 static void createFieldNames
00225 (
00226     const Time& runTime,
00227     const instantList& Times,
00228     const word& setName
00229 )
00230 {
00231     // From foamToFieldView9/getFieldNames.H:
00232 
00233     HashSet<word> volScalarHash;
00234     HashSet<word> volVectorHash;
00235     HashSet<word> surfScalarHash;
00236     HashSet<word> surfVectorHash;
00237 
00238     if (setName.empty())
00239     {
00240         forAll(Times, timeI)
00241         {
00242             const word& timeName = Times[timeI].name();
00243 
00244             // Add all fields to hashtable
00245             IOobjectList objects(runTime, timeName);
00246 
00247             wordList vsNames(objects.names(volScalarField::typeName));
00248 
00249             forAll(vsNames, fieldI)
00250             {
00251                 volScalarHash.insert(vsNames[fieldI]);
00252             }
00253 
00254             wordList vvNames(objects.names(volVectorField::typeName));
00255 
00256             forAll(vvNames, fieldI)
00257             {
00258                 volVectorHash.insert(vvNames[fieldI]);
00259             }
00260         }
00261     }
00262     db_.setFieldNames(volScalarHash.toc(), volVectorHash.toc());
00263 }
00264 
00265 
00266 // Appends interpolated values of fieldName to vars array.
00267 static void storeScalarField
00268 (
00269     const volPointInterpolation& pInterp,
00270     const word& fieldName,
00271     float vars[],
00272     label& pointI
00273 )
00274 {
00275     label nPoints = db_.mesh().nPoints();
00276     label nTotPoints = nPoints + db_.polys().size();
00277 
00278     // Check if present
00279     IOobject ioHeader
00280     (
00281         fieldName,
00282         db_.runTime().timeName(),
00283         db_.runTime(),
00284         IOobject::MUST_READ,
00285         IOobject::NO_WRITE
00286     );
00287 
00288     if (ioHeader.headerOk())
00289     {
00290         Info<< "Storing " << nTotPoints << " of interpolated " << fieldName
00291             << endl;
00292 
00293         volScalarField field(ioHeader, db_.mesh());
00294 
00295         pointScalarField psf(pInterp.interpolate(field));
00296 
00297         forAll(psf, i)
00298         {
00299             vars[pointI++] = float(psf[i]);
00300         }
00301 
00302         const labelList& polys = db_.polys();
00303 
00304         forAll(polys, i)
00305         {
00306             label cellI = polys[i];
00307 
00308             vars[pointI++] = float(field[cellI]);
00309         }
00310     }
00311     else
00312     {
00313         Info<< "Storing " << nTotPoints << " of dummy values of " << fieldName
00314             << endl;
00315 
00316         for(label i = 0; i < nPoints; i++)
00317         {
00318             vars[pointI++] = 0.0;
00319         }
00320 
00321         const labelList& polys = db_.polys();
00322 
00323         forAll(polys, i)
00324         {
00325             vars[pointI++] = 0.0;
00326         }
00327     }
00328 }
00329 
00330 
00331 // Appends interpolated values of fieldName to vars array.
00332 static void storeVectorField
00333 (
00334     const volPointInterpolation& pInterp,
00335     const word& fieldName,
00336     float vars[],
00337     label& pointI
00338 )
00339 {
00340     label nPoints = db_.mesh().nPoints();
00341 
00342     label nTotPoints = nPoints + db_.polys().size();
00343 
00344     // Check if present
00345     IOobject ioHeader
00346     (
00347         fieldName,
00348         db_.runTime().timeName(),
00349         db_.runTime(),
00350         IOobject::MUST_READ,
00351         IOobject::NO_WRITE
00352     );
00353 
00354     if (ioHeader.headerOk())
00355     {
00356         Info<< "Storing " << nTotPoints << " of interpolated " << fieldName
00357             << endl;
00358 
00359         volVectorField field(ioHeader, db_.mesh());
00360 
00361         for (direction d = 0; d < vector::nComponents; d++)
00362         {
00363             tmp<volScalarField> tcomp = field.component(d);
00364             const volScalarField& comp = tcomp();
00365 
00366             pointScalarField psf(pInterp.interpolate(comp));
00367 
00368             forAll(psf, i)
00369             {
00370                 vars[pointI++] = float(psf[i]);
00371             }
00372 
00373             const labelList& polys = db_.polys();
00374 
00375             forAll(polys, i)
00376             {
00377                 label cellI = polys[i];
00378 
00379                 vars[pointI++] = float(comp[cellI]);
00380             }
00381         }
00382     }
00383     else
00384     {
00385         Info<< "Storing " << nTotPoints << " of dummy values of " << fieldName
00386             << endl;
00387 
00388         for (direction d = 0; d < vector::nComponents; d++)
00389         {
00390             for(label i = 0; i < nPoints; i++)
00391             {
00392                 vars[pointI++] = 0.0;
00393             }
00394 
00395             const labelList& polys = db_.polys();
00396 
00397             forAll(polys, i)
00398             {
00399                 vars[pointI++] = 0.0;
00400             }
00401         }
00402     }
00403 }
00404 
00405 
00406 // Returns Fieldview face_type of mesh face faceI.
00407 static label getFvType(const polyMesh& mesh, const label faceI)
00408 {
00409     return mesh.boundaryMesh().whichPatch(faceI) + 1;
00410 }
00411 
00412 
00413 // Returns Fieldview face_type of face f.
00414 static label getFaceType
00415 (
00416     const polyMesh& mesh,
00417     const labelList& faceLabels,
00418     const face& f
00419 )
00420 {
00421     // Search in subset faceLabels of faces for index of face f.
00422     const faceList& faces = mesh.faces();
00423 
00424     forAll(faceLabels, i)
00425     {
00426         label faceI = faceLabels[i];
00427 
00428         if (f == faces[faceI])
00429         {
00430             // Convert patch to Fieldview face_type.
00431             return getFvType(mesh, faceI);
00432         }
00433     }
00434 
00435     FatalErrorIn("getFaceType")
00436         << "Cannot find face " << f << " in mesh face subset " << faceLabels
00437         << abort(FatalError);
00438 
00439     return -1;
00440 }
00441 
00442 
00443 // Returns Fieldview face_types for set of faces
00444 static labelList getFaceTypes
00445 (
00446     const polyMesh& mesh,
00447     const labelList& cellFaces,
00448     const faceList& cellShapeFaces
00449 )
00450 {
00451     labelList faceLabels(cellShapeFaces.size());
00452 
00453     forAll(cellShapeFaces, i)
00454     {
00455         faceLabels[i] = getFaceType(mesh, cellFaces, cellShapeFaces[i]);
00456     }
00457     return faceLabels;
00458 }
00459 
00460 
00461 /*
00462  * Callback for querying file contents. Taken from user/user_unstruct_combined.f
00463  */
00464 void user_query_file_function
00465 (
00466     /* input */
00467     char* fname,        /* filename */
00468     int* lenf,          /* length of fName */
00469     int* iunit,         /* fortran unit to use */
00470     int* max_grids,     /* maximum number of grids allowed */
00471     int* max_face_types,/* maximum number of face types allowed */
00472     int* max_vars,      /* maximum number of result variables allowed per */
00473                         /* grid point*/
00474 
00475     /* output */
00476     int* num_grids,     /* number of grids that will be read from data file */
00477     int  num_nodes[],   /* (array of node counts for all grids) */
00478     int* num_face_types,        /* number of special face types */
00479     char face_type_names[][80], /* array of face-type names */
00480     int  wall_flags[],          /* array of flags for the "wall" behavior */
00481     int* num_vars,              /* number of result variables per grid point */
00482     char var_names[][80],       /* array of variable names */
00483     int* iret                   /* return value */
00484 )
00485 {
00486     fprintf(stderr, "\n** user_query_file_function\n");
00487 
00488     string rootAndCaseString(fname, *lenf);
00489 
00490     fileName rootAndCase(rootAndCaseString);
00491 
00492     word setName("");
00493 
00494     if (!validCase(rootAndCase))
00495     {
00496         setName = rootAndCase.name();
00497 
00498         rootAndCase = rootAndCase.path();
00499 
00500         word setDir = rootAndCase.name();
00501 
00502         rootAndCase = rootAndCase.path();
00503 
00504         word meshDir = rootAndCase.name();
00505 
00506         rootAndCase = rootAndCase.path();
00507         rootAndCase = rootAndCase.path();
00508 
00509         if
00510         (
00511             setDir == "sets"
00512          && meshDir == polyMesh::typeName
00513          && validCase(rootAndCase)
00514         )
00515         {
00516             // Valid set (hopefully - cannot check contents of setName yet).
00517         }
00518         else
00519         {
00520             errorMsg
00521             (
00522                 "Could not find system/ and constant/ directory in\n"
00523               + rootAndCase
00524               + "\nPlease select a Foam case directory."
00525             );
00526 
00527             *iret = 1;
00528 
00529             return;
00530         }
00531 
00532     }
00533 
00534     fileName rootDir(rootAndCase.path());
00535 
00536     fileName caseName(rootAndCase.name());
00537 
00538     // handle trailing '/'
00539     if (caseName.empty())
00540     {
00541         caseName = rootDir.name();
00542         rootDir  = rootDir.path();
00543     }
00544 
00545     Info<< "rootDir  : " << rootDir << endl
00546         << "caseName : " << caseName << endl
00547         << "setName  : " << setName << endl;
00548 
00549     //
00550     // Get/reuse database and mesh
00551     //
00552 
00553     bool caseChanged = db_.setRunTime(rootDir, caseName, setName);
00554 
00555 
00556     //
00557     // Select time
00558     //
00559 
00560     instantList Times = db_.runTime().times();
00561 
00562     // If topo case set database time and update mesh.
00563     if (hasTopoChange(Times))
00564     {
00565         if (!selectTime(Times, iret))
00566         {
00567             return;
00568         }
00569     }
00570     else if (caseChanged)
00571     {
00572         // Load mesh (if case changed) to make sure we have nPoints etc.
00573         db_.loadMesh();
00574     }
00575 
00576 
00577     //
00578     // Set output variables
00579     //
00580 
00581     *num_grids = 1;
00582 
00583     const fvMesh& mesh = db_.mesh();
00584 
00585     label nTotPoints = mesh.nPoints() + db_.polys().size();
00586 
00587     num_nodes[0] = nTotPoints;
00588 
00589     Info<< "setting num_nodes:" << num_nodes[0] << endl;
00590 
00591     Info<< "setting num_face_types:" << mesh.boundary().size() << endl;
00592 
00593     *num_face_types = mesh.boundary().size();
00594 
00595     if (*num_face_types > *max_face_types)
00596     {
00597         errorMsg("Too many patches. FieldView limit:" + name(*max_face_types));
00598 
00599         *iret = 1;
00600 
00601         return;
00602     }
00603 
00604 
00605     forAll(mesh.boundaryMesh(), patchI)
00606     {
00607         const polyPatch& patch = mesh.boundaryMesh()[patchI];
00608 
00609         strcpy(face_type_names[patchI], patch.name().c_str());
00610 
00611         if (isA<wallPolyPatch>(patch))
00612         {
00613             wall_flags[patchI] = 1;
00614         }
00615         else
00616         {
00617             wall_flags[patchI] = 0;
00618         }
00619         Info<< "Patch " << patch.name() << " is wall:"
00620             <<  wall_flags[patchI] << endl;
00621     }
00622 
00623     //- Find all volFields and add them to database
00624     createFieldNames(db_.runTime(), Times, setName);
00625 
00626     *num_vars = db_.volScalarNames().size() + 3*db_.volVectorNames().size();
00627 
00628     if (*num_vars > *max_vars)
00629     {
00630         errorMsg("Too many variables. FieldView limit:" + name(*max_vars));
00631 
00632         *iret = 1;
00633 
00634         return;
00635     }
00636 
00637 
00638     int nameI = 0;
00639 
00640     forAll(db_.volScalarNames(), i)
00641     {
00642         const word& fieldName = db_.volScalarNames()[i];
00643 
00644         const word& fvName = db_.getFvName(fieldName);
00645 
00646         strcpy(var_names[nameI++], fvName.c_str());
00647     }
00648 
00649     forAll(db_.volVectorNames(), i)
00650     {
00651         const word& fieldName = db_.volVectorNames()[i];
00652 
00653         const word& fvName = db_.getFvName(fieldName);
00654 
00655         strcpy(var_names[nameI++], (fvName + "x;" + fvName).c_str());
00656         strcpy(var_names[nameI++], (fvName + "y").c_str());
00657         strcpy(var_names[nameI++], (fvName + "z").c_str());
00658     }
00659 
00660     *iret = 0;
00661 }
00662 
00663 
00664 /*
00665  * Callback for reading timestep. Taken from user/user_unstruct_combined.f
00666  */
00667 void user_read_one_grid_function
00668 (
00669     int* iunit,     /* in: fortran unit number */
00670     int* igrid,     /* in: grid number to read */
00671     int* nodecnt,   /* in: number of nodes to read */
00672     float xyz[],    /* out: coordinates of nodes: x1..xN y1..yN z1..zN */
00673     int* num_vars,  /* in: number of results per node */
00674     float vars[],   /* out: values per node */
00675     int* iret       /* out: return value */
00676 )
00677 {
00678     fprintf(stderr, "\n** user_read_one_grid_function\n");
00679 
00680     if (*igrid != 1)
00681     {
00682         errorMsg("Illegal grid number " + Foam::name(*igrid));
00683 
00684         *iret = 1;
00685 
00686         return;
00687     }
00688 
00689     // Get current time
00690     instantList Times = db_.runTime().times();
00691 
00692     // Set database time and update mesh.
00693     // Note: this should not be nessecary here. We already have the correct
00694     // time set and mesh loaded. This is only nessecary because Fieldview
00695     // otherwise thinks the case is non-transient.
00696     if (!selectTime(Times, iret))
00697     {
00698         return;
00699     }
00700 
00701 
00702     const fvMesh& mesh = db_.mesh();
00703 
00704     // With mesh now loaded check for change in number of points.
00705     label nTotPoints = mesh.nPoints() + db_.polys().size();
00706 
00707     if (*nodecnt != nTotPoints)
00708     {
00709         errorMsg
00710         (
00711             "nodecnt differs from number of points in mesh.\nnodecnt:"
00712           + Foam::name(*nodecnt)
00713           + "  mesh:"
00714           + Foam::name(nTotPoints)
00715         );
00716 
00717         *iret = 1;
00718 
00719         return;
00720     }
00721 
00722 
00723     if
00724     (
00725         *num_vars
00726      != (db_.volScalarNames().size() + 3*db_.volVectorNames().size())
00727     )
00728     {
00729         errorMsg("Illegal number of variables " + name(*num_vars));
00730 
00731         *iret = 1;
00732 
00733         return;
00734     }
00735 
00736     //
00737     // Set coordinates
00738     //
00739 
00740     const pointField& points = mesh.points();
00741 
00742     int xIndex = 0;
00743     int yIndex = xIndex + nTotPoints;
00744     int zIndex = yIndex + nTotPoints;
00745 
00746     // Add mesh points first.
00747     forAll(points, pointI)
00748     {
00749         xyz[xIndex++] = points[pointI].x();
00750         xyz[yIndex++] = points[pointI].y();
00751         xyz[zIndex++] = points[pointI].z();
00752     }
00753 
00754     // Add cell centres of polys
00755     const pointField& ctrs = mesh.cellCentres();
00756 
00757     const labelList& polys = db_.polys();
00758 
00759     forAll(polys, i)
00760     {
00761         label cellI = polys[i];
00762 
00763         xyz[xIndex++] = ctrs[cellI].x();
00764         xyz[yIndex++] = ctrs[cellI].y();
00765         xyz[zIndex++] = ctrs[cellI].z();
00766     }
00767 
00768 
00769     //
00770     // Define elements by calling fv routines
00771     //
00772 
00773     static const cellModel& tet = *(cellModeller::lookup("tet"));
00774     static const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
00775     static const cellModel& pyr = *(cellModeller::lookup("pyr"));
00776     static const cellModel& prism = *(cellModeller::lookup("prism"));
00777     static const cellModel& wedge = *(cellModeller::lookup("wedge"));
00778     static const cellModel& hex = *(cellModeller::lookup("hex"));
00779     //static const cellModel& splitHex = *(cellModeller::lookup("splitHex"));
00780 
00781     int tetVerts[4];
00782     int pyrVerts[5];
00783     int prismVerts[6];
00784     int hexVerts[8];
00785 
00786     int tetFaces[4];
00787     int pyrFaces[5];
00788     int prismFaces[5];
00789     int hexFaces[6];
00790 
00791     const cellShapeList& cellShapes = mesh.cellShapes();
00792     const faceList& faces = mesh.faces();
00793     const cellList& cells = mesh.cells();
00794     const labelList& owner = mesh.faceOwner();
00795     label nPoints = mesh.nPoints();
00796 
00797     // Fieldview face_types array with all faces marked as internal.
00798     labelList internalFaces(6, 0);
00799 
00800     // Mark all cells next to boundary so we don't have to calculate
00801     // wall_types for internal cells and can just pass internalFaces.
00802     boolList wallCell(mesh.nCells(), false);
00803 
00804     label nWallCells = 0;
00805 
00806     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
00807     {
00808         label cellI = owner[faceI];
00809 
00810         if (!wallCell[cellI])
00811         {
00812             wallCell[cellI] = true;
00813 
00814             nWallCells++;
00815         }
00816     }
00817 
00818     label nPolys = 0;
00819 
00820     forAll(cellShapes, cellI)
00821     {
00822         const cellShape& cellShape = cellShapes[cellI];
00823         const cellModel& cellModel = cellShape.model();
00824         const cell& cellFaces = cells[cellI];
00825 
00826         int istat = 0;
00827 
00828         if (cellModel == tet)
00829         {
00830             tetVerts[0] = cellShape[3] + 1;
00831             tetVerts[1] = cellShape[0] + 1;
00832             tetVerts[2] = cellShape[1] + 1;
00833             tetVerts[3] = cellShape[2] + 1;
00834 
00835             if (wallCell[cellI])
00836             {
00837                 labelList faceTypes =
00838                     getFaceTypes(mesh, cellFaces, cellShape.faces());
00839 
00840                 tetFaces[0] = faceTypes[2];
00841                 tetFaces[1] = faceTypes[3];
00842                 tetFaces[2] = faceTypes[0];
00843                 tetFaces[3] = faceTypes[1];
00844 
00845                 istat = create_tet(ITYPE, tetVerts, tetFaces);
00846             }
00847             else
00848             {
00849                 // All faces internal so use precalculated zero.
00850                 istat = create_tet(ITYPE, tetVerts, internalFaces.begin());
00851             }
00852         }
00853         else if (cellModel == tetWedge)
00854         {
00855             prismVerts[0] = cellShape[0] + 1;
00856             prismVerts[1] = cellShape[3] + 1;
00857             prismVerts[2] = cellShape[4] + 1;
00858             prismVerts[3] = cellShape[1] + 1;
00859             prismVerts[4] = cellShape[4] + 1;
00860             prismVerts[5] = cellShape[2] + 1;
00861 
00862             if (wallCell[cellI])
00863             {
00864                 labelList faceTypes =
00865                     getFaceTypes(mesh, cellFaces, cellShape.faces());
00866 
00867                 prismFaces[0] = faceTypes[1];
00868                 prismFaces[1] = faceTypes[2];
00869                 prismFaces[2] = faceTypes[3];
00870                 prismFaces[3] = faceTypes[0];
00871                 prismFaces[4] = faceTypes[3];
00872 
00873                 istat = create_prism(ITYPE, prismVerts, prismFaces);
00874             }
00875             else
00876             {
00877                 istat = create_prism(ITYPE, prismVerts, internalFaces.begin());
00878             }
00879         }
00880         else if (cellModel == pyr)
00881         {
00882             pyrVerts[0] = cellShape[0] + 1;
00883             pyrVerts[1] = cellShape[1] + 1;
00884             pyrVerts[2] = cellShape[2] + 1;
00885             pyrVerts[3] = cellShape[3] + 1;
00886             pyrVerts[4] = cellShape[4] + 1;
00887 
00888             if (wallCell[cellI])
00889             {
00890                 labelList faceTypes =
00891                     getFaceTypes(mesh, cellFaces, cellShape.faces());
00892 
00893                 pyrFaces[0] = faceTypes[0];
00894                 pyrFaces[1] = faceTypes[3];
00895                 pyrFaces[2] = faceTypes[2];
00896                 pyrFaces[3] = faceTypes[1];
00897                 pyrFaces[4] = faceTypes[4];
00898 
00899                 istat = create_pyramid(ITYPE, pyrVerts, pyrFaces);
00900             }
00901             else
00902             {
00903                 istat = create_pyramid(ITYPE, pyrVerts, internalFaces.begin());
00904             }
00905         }
00906         else if (cellModel == prism)
00907         {
00908             prismVerts[0] = cellShape[0] + 1;
00909             prismVerts[1] = cellShape[3] + 1;
00910             prismVerts[2] = cellShape[4] + 1;
00911             prismVerts[3] = cellShape[1] + 1;
00912             prismVerts[4] = cellShape[5] + 1;
00913             prismVerts[5] = cellShape[2] + 1;
00914 
00915             if (wallCell[cellI])
00916             {
00917                 labelList faceTypes =
00918                     getFaceTypes(mesh, cellFaces, cellShape.faces());
00919 
00920                 prismFaces[0] = faceTypes[4];
00921                 prismFaces[1] = faceTypes[2];
00922                 prismFaces[2] = faceTypes[3];
00923                 prismFaces[3] = faceTypes[0];
00924                 prismFaces[4] = faceTypes[1];
00925 
00926                 istat = create_prism(ITYPE, prismVerts, prismFaces);
00927             }
00928             else
00929             {
00930                 istat = create_prism(ITYPE, prismVerts, internalFaces.begin());
00931             }
00932         }
00933         else if (cellModel == wedge)
00934         {
00935             hexVerts[0] = cellShape[0] + 1;
00936             hexVerts[1] = cellShape[1] + 1;
00937             hexVerts[2] = cellShape[0] + 1;
00938             hexVerts[3] = cellShape[2] + 1;
00939             hexVerts[4] = cellShape[3] + 1;
00940             hexVerts[5] = cellShape[4] + 1;
00941             hexVerts[6] = cellShape[6] + 1;
00942             hexVerts[7] = cellShape[5] + 1;
00943 
00944             if (wallCell[cellI])
00945             {
00946                 labelList faceTypes =
00947                     getFaceTypes(mesh, cellFaces, cellShape.faces());
00948 
00949                 hexFaces[0] = faceTypes[2];
00950                 hexFaces[1] = faceTypes[3];
00951                 hexFaces[2] = faceTypes[0];
00952                 hexFaces[3] = faceTypes[1];
00953                 hexFaces[4] = faceTypes[4];
00954                 hexFaces[5] = faceTypes[5];
00955 
00956                 istat = create_hex(ITYPE, hexVerts, hexFaces);
00957             }
00958             else
00959             {
00960                 istat = create_hex(ITYPE, hexVerts, internalFaces.begin());
00961             }
00962         }
00963         else if (cellModel == hex)
00964         {
00965             hexVerts[0] = cellShape[0] + 1;
00966             hexVerts[1] = cellShape[1] + 1;
00967             hexVerts[2] = cellShape[3] + 1;
00968             hexVerts[3] = cellShape[2] + 1;
00969             hexVerts[4] = cellShape[4] + 1;
00970             hexVerts[5] = cellShape[5] + 1;
00971             hexVerts[6] = cellShape[7] + 1;
00972             hexVerts[7] = cellShape[6] + 1;
00973 
00974             if (wallCell[cellI])
00975             {
00976                 labelList faceTypes =
00977                     getFaceTypes(mesh, cellFaces, cellShape.faces());
00978 
00979                 hexFaces[0] = faceTypes[0];
00980                 hexFaces[1] = faceTypes[1];
00981                 hexFaces[2] = faceTypes[4];
00982                 hexFaces[3] = faceTypes[5];
00983                 hexFaces[4] = faceTypes[2];
00984                 hexFaces[5] = faceTypes[3];
00985 
00986                 istat = create_hex(ITYPE, hexVerts, hexFaces);
00987             }
00988             else
00989             {
00990                 istat = create_hex(ITYPE, hexVerts, internalFaces.begin());
00991             }
00992         }
00993         else
00994         {
00995             forAll(cellFaces, cFaceI)
00996             {
00997                 label faceI = cellFaces[cFaceI];
00998 
00999                 // Get Fieldview facetype (internal/on patch)
01000                 label fvFaceType = getFvType(mesh, faceI);
01001 
01002                 const face& f = faces[faceI];
01003 
01004                 // Indices into storage
01005                 label nQuads = 0;
01006                 label nTris = 0;
01007 
01008                 // Storage for triangles and quads created by face
01009                 // decomposition (sized for worst case)
01010                 faceList quadFaces((f.size() - 2)/2);
01011                 faceList triFaces(f.size() - 2);
01012 
01013                 f.trianglesQuads
01014                 (
01015                     points,
01016                     nTris,
01017                     nQuads,
01018                     triFaces,
01019                     quadFaces
01020                 );
01021 
01022                 // Label of cell centre in fv point list.
01023                 label polyCentrePoint = nPoints + nPolys;
01024 
01025                 for (label i=0; i<nTris; i++)
01026                 {
01027                     if (cellI == owner[faceI])
01028                     {
01029                         tetVerts[0] = triFaces[i][0] + 1;
01030                         tetVerts[1] = triFaces[i][1] + 1;
01031                         tetVerts[2] = triFaces[i][2] + 1;
01032                         tetVerts[3] = polyCentrePoint + 1;
01033                     }
01034                     else
01035                     {
01036                         tetVerts[0] = triFaces[i][2] + 1;
01037                         tetVerts[1] = triFaces[i][1] + 1;
01038                         tetVerts[2] = triFaces[i][0] + 1;
01039                         tetVerts[3] = polyCentrePoint + 1;
01040                     }
01041 
01042                     if (wallCell[cellI])
01043                     {
01044                         // Outside face is one without polyCentrePoint
01045                         tetFaces[0] = fvFaceType;
01046                         tetFaces[1] = 0;
01047                         tetFaces[2] = 0;
01048                         tetFaces[3] = 0;
01049 
01050                         istat = create_tet(ITYPE, tetVerts, tetFaces);
01051                     }
01052                     else
01053                     {
01054                         istat =
01055                             create_tet
01056                             (
01057                                 ITYPE,
01058                                 tetVerts,
01059                                 internalFaces.begin()
01060                             );
01061                     }
01062                 }
01063 
01064                 for (label i=0; i<nQuads; i++)
01065                 {
01066                     if (cellI == owner[faceI])
01067                     {
01068                         pyrVerts[0] = quadFaces[i][3] + 1;
01069                         pyrVerts[1] = quadFaces[i][2] + 1;
01070                         pyrVerts[2] = quadFaces[i][1] + 1;
01071                         pyrVerts[3] = quadFaces[i][0] + 1;
01072                         pyrVerts[4] = polyCentrePoint + 1;
01073 
01074                     }
01075                     else
01076                     {
01077                         pyrVerts[0] = quadFaces[i][0] + 1;
01078                         pyrVerts[1] = quadFaces[i][1] + 1;
01079                         pyrVerts[2] = quadFaces[i][2] + 1;
01080                         pyrVerts[3] = quadFaces[i][3] + 1;
01081                         pyrVerts[4] = polyCentrePoint + 1;
01082                     }
01083 
01084                     if (wallCell[cellI])
01085                     {
01086                         // Outside face is one without polyCentrePoint
01087                         pyrFaces[0] = fvFaceType;
01088                         pyrFaces[1] = 0;
01089                         pyrFaces[2] = 0;
01090                         pyrFaces[3] = 0;
01091                         pyrFaces[4] = 0;
01092 
01093                         istat = create_pyramid(ITYPE, pyrVerts, pyrFaces);
01094                     }
01095                     else
01096                     {
01097                         istat =
01098                             create_pyramid
01099                             (
01100                                 ITYPE,
01101                                 pyrVerts,
01102                                 internalFaces.begin()
01103                             );
01104                     }
01105                 }
01106 
01107                 if (istat != 0)
01108                 {
01109                     errorMsg("Error during adding cell " + name(cellI));
01110 
01111                     *iret = 1;
01112 
01113                     return;
01114                 }
01115             }
01116 
01117             nPolys++;
01118         }
01119 
01120         if (istat != 0)
01121         {
01122             errorMsg("Error during adding cell " + name(cellI));
01123 
01124             *iret = 1;
01125 
01126             return;
01127         }
01128     }
01129 
01130 
01131 
01132     //
01133     // Set fieldvalues
01134     //
01135 
01136     pointMesh pMesh(mesh);
01137 
01138     volPointInterpolation pInterp(mesh, pMesh);
01139 
01140     int pointI = 0;
01141 
01142     forAll(db_.volScalarNames(), i)
01143     {
01144         const word& fieldName = db_.volScalarNames()[i];
01145 
01146         storeScalarField(pInterp, fieldName, vars, pointI);
01147     }
01148 
01149 
01150     forAll(db_.volVectorNames(), i)
01151     {
01152         const word& fieldName = db_.volVectorNames()[i];
01153 
01154         storeVectorField(pInterp, fieldName, vars, pointI);
01155     }
01156 
01157     // Return without failure.
01158     *iret = 0;
01159 }
01160 
01161 
01162 void register_data_readers()
01163 {
01164     /*
01165     **
01166     ** You should edit this file to "register" a user-defined data
01167     ** reader with FIELDVIEW, if the user functions being registered
01168     ** here are written in C.
01169     ** You should edit "ftn_register_data_readers.f" if the user functions
01170     ** being registered are written in Fortran.
01171     ** In either case, the user functions being registered may call other
01172     ** functions written in either language (C or Fortran); only the
01173     ** language of the "query" and "read" functions referenced here matters
01174     ** to FIELDVIEW.
01175     **
01176     ** The following shows a sample user-defined data reader being
01177     ** "registered" with FIELDVIEW.
01178     **
01179     ** The "extern void" declarations should match the names of the
01180     ** query and grid-reading functions you are providing. It is
01181     ** strongly suggested that all such names begin with "user" so
01182     ** as not to conflict with global names in FIELDVIEW.
01183     **
01184     ** You may call any combination of the data reader registration
01185     ** functions shown below ("register_two_file_reader" and/or
01186     ** "register_single_file_reader" and/or "register_single_unstruct_reader"
01187     ** and/or "register_double_unstruct_reader") as many times as you like,
01188     ** in order to create several different data readers. Each data reader
01189     ** should of course have different "query" and "read" functions, all of
01190     ** which should also appear in "extern void" declarations.
01191     **
01192     */
01193 
01194     /*
01195     ** like this for combined unstructured grids & results in a single file
01196     */
01197     reg_single_unstruct_reader (
01198     "Foam Reader",             /* title you want for data reader */
01199     user_query_file_function,      /* whatever you called this */
01200     user_read_one_grid_function    /* whatever you called this */
01201     );
01202 }
01203 
01204 
01205 
01206 
01207 }
01208 
01209 
01210 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines