avtExtrudedVolFileFormat.C 11.9 KB
Newer Older
1 2
/*****************************************************************************
*
bonnell's avatar
bonnell committed
3
* Copyright (c) 2000 - 2015, Lawrence Livermore National Security, LLC
4
* Produced at the Lawrence Livermore National Laboratory
5
* LLNL-CODE-442911
6 7
* All rights reserved.
*
8
* This file is  part of VisIt. For  details, see https://visit.llnl.gov/.  The
9 10 11 12 13 14 15 16 17 18
* full copyright notice is contained in the file COPYRIGHT located at the root
* of the VisIt distribution or at http://www.llnl.gov/visit/copyright.html.
*
* Redistribution  and  use  in  source  and  binary  forms,  with  or  without
* modification, are permitted provided that the following conditions are met:
*
*  - Redistributions of  source code must  retain the above  copyright notice,
*    this list of conditions and the disclaimer below.
*  - Redistributions in binary form must reproduce the above copyright notice,
*    this  list of  conditions  and  the  disclaimer (as noted below)  in  the
19 20 21
*    documentation and/or other materials provided with the distribution.
*  - Neither the name of  the LLNS/LLNL nor the names of  its contributors may
*    be used to endorse or promote products derived from this software without
22 23 24 25 26
*    specific prior written permission.
*
* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT  HOLDERS AND CONTRIBUTORS "AS IS"
* AND ANY EXPRESS OR  IMPLIED WARRANTIES, INCLUDING,  BUT NOT  LIMITED TO, THE
* IMPLIED WARRANTIES OF MERCHANTABILITY AND  FITNESS FOR A PARTICULAR  PURPOSE
27 28 29
* ARE  DISCLAIMED. IN  NO EVENT  SHALL LAWRENCE  LIVERMORE NATIONAL  SECURITY,
* LLC, THE  U.S.  DEPARTMENT OF  ENERGY  OR  CONTRIBUTORS BE  LIABLE  FOR  ANY
* DIRECT,  INDIRECT,   INCIDENTAL,   SPECIAL,   EXEMPLARY,  OR   CONSEQUENTIAL
30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45
* DAMAGES (INCLUDING, BUT NOT  LIMITED TO, PROCUREMENT OF  SUBSTITUTE GOODS OR
* SERVICES; LOSS OF  USE, DATA, OR PROFITS; OR  BUSINESS INTERRUPTION) HOWEVER
* CAUSED  AND  ON  ANY  THEORY  OF  LIABILITY,  WHETHER  IN  CONTRACT,  STRICT
* LIABILITY, OR TORT  (INCLUDING NEGLIGENCE OR OTHERWISE)  ARISING IN ANY  WAY
* OUT OF THE  USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH
* DAMAGE.
*
*****************************************************************************/

// ************************************************************************* //
//                         avtExtrudedVolFileFormat.C                        //
// ************************************************************************* //

#include <avtExtrudedVolFileFormat.h>

#include <string>
bonnell's avatar
bonnell committed
46
#include <visitstream.h>
47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169
#include <snprintf.h>

#include <vtkCell.h>
#include <vtkCellTypes.h>
#include <vtkDataSetReader.h>
#include <vtkFloatArray.h>
#include <vtkRectilinearGrid.h>
#include <vtkUnstructuredGrid.h>

#include <avtDatabaseMetaData.h>

#include <Expression.h>

#include <InvalidVariableException.h>
#include <InvalidFilesException.h>


using     std::string;


// ****************************************************************************
//  Method: avtExtrudedVol constructor
//
//  Programmer: childs -- generated by xml2avt
//  Creation:   Fri May 18 17:52:04 PST 2007
//
// ****************************************************************************

avtExtrudedVolFileFormat::avtExtrudedVolFileFormat(const char *filename, DBOptionsAttributes *readOpts)
    : avtSTMDFileFormat(&filename, 1)
{
    ifstream ifile(filename);
    if (ifile.fail())
    {
        EXCEPTION1(InvalidFilesException, filename);
    }
 
    std::string tmp;
    ifile >> tmp; // "STEM: ";
    ifile >> tmp;
    stem = tmp;

    ifile >> tmp; // "NUMCHUNKS: ";
    ifile >> tmp;
    numChunks = atoi(tmp.c_str());
    if (numChunks <= 0)
    {
        EXCEPTION1(VisItException, "This does not appear to be a valid "
                                   "ExtrudedVol file.");
    }

    ifile >> tmp; // "NTIMES: ";
    ifile >> tmp;
    nTimesteps = atoi(tmp.c_str());
    if (nTimesteps <= 0)
    {
        EXCEPTION1(VisItException, "This does not appear to be a valid "
                                   "ExtrudedVol file.");
    }

    ifile >> tmp; // "VARIABLES: ";
    ifile >> tmp;
    int nVars = atoi(tmp.c_str());
    if (nVars <= 0)
    {
        EXCEPTION1(VisItException, "This does not appear to be a valid "
                                   "ExtrudedVol file.");
    }
    for (int i = 0 ; i < nVars ; i++)
    {
        ifile >> tmp;
        variables.push_back(tmp);
    }
}


// ****************************************************************************
//  Method: avtExtrudedVolFileFormat::FreeUpResources
//
//  Purpose:
//      When VisIt is done focusing on a particular timestep, it asks that
//      timestep to free up any resources (memory, file descriptors) that
//      it has associated with it.  This method is the mechanism for doing
//      that.
//
//  Programmer: childs -- generated by xml2avt
//  Creation:   Fri May 18 17:52:04 PST 2007
//
// ****************************************************************************

void
avtExtrudedVolFileFormat::FreeUpResources(void)
{
}


// ****************************************************************************
//  Method: avtExtrudedVolFileFormat::PopulateDatabaseMetaData
//
//  Purpose:
//      This database meta-data object is like a table of contents for the
//      file.  By populating it, you are telling the rest of VisIt what
//      information it can request from you.
//
//  Programmer: childs -- generated by xml2avt
//  Creation:   Fri May 18 17:52:04 PST 2007
//
// ****************************************************************************

void
avtExtrudedVolFileFormat::PopulateDatabaseMetaData(avtDatabaseMetaData *md)
{
    string meshname = "wedge_mesh";
    avtMeshType mt = AVT_UNSTRUCTURED_MESH;
    int nblocks = numChunks;
    int block_origin = 0;
    int spatial_dimension = 3;
    int topological_dimension = 3;
    double *extents = NULL;
    AddMeshToMetaData(md, meshname, mt, extents, nblocks, block_origin,
                      spatial_dimension, topological_dimension);

    string mesh_for_this_var = meshname; 
170
    for (size_t i = 0 ; i < variables.size() ; i++)
171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204
        AddScalarVarToMetaData(md, variables[i], mesh_for_this_var, 
                               AVT_NODECENT);
}


// ****************************************************************************
//  Method: avtExtrudedVolFileFormat::GetMesh
//
//  Purpose:
//      Gets the mesh associated with this file.  The mesh is returned as a
//      derived type of vtkDataSet (ie vtkRectilinearGrid, vtkStructuredGrid,
//      vtkUnstructuredGrid, etc).
//
//  Arguments:
//      domain      The index of the domain.  If there are NDomains, this
//                  value is guaranteed to be between 0 and NDomains-1,
//                  regardless of block origin.
//      meshname    The name of the mesh of interest.  This can be ignored if
//                  there is only one mesh.
//
//  Programmer: childs -- generated by xml2avt
//  Creation:   Fri May 18 17:52:04 PST 2007
//
// ****************************************************************************

vtkDataSet *
avtExtrudedVolFileFormat::GetMesh(int domain, const char *meshname)
{
    int   i, j;

    char filename[1024];
    SNPRINTF(filename, 1024, "%s.%d.exvol_conn", stem.c_str(), domain);
    vtkDataSetReader *rdr = vtkDataSetReader::New();
    rdr->SetFileName(filename);
bonnell's avatar
bonnell committed
205
    rdr->Update();
206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247
    vtkDataSet *ds = rdr->GetOutput();
    if (ds == NULL || ds->GetDataObjectType() != VTK_POLY_DATA)
    {
        EXCEPTION1(InvalidVariableException, meshname);
    }

    vtkUnstructuredGrid *ugrid = vtkUnstructuredGrid::New();
    vtkPoints *pts = vtkPoints::New();
    int npts     = ds->GetNumberOfPoints()*nTimesteps;
    int nOrigPts = ds->GetNumberOfPoints();
    pts->SetNumberOfPoints(npts);

    for (i = 0 ; i < nTimesteps ; i++)
    {
        SNPRINTF(filename, 1024, "%s.%d.%d.exvol_var", stem.c_str(), i,domain);
        ifstream ifile(filename);
        if (ifile.fail())
        {
            EXCEPTION1(InvalidVariableException, meshname);
        }
        int npts;
        ifile >> npts;
        for (j = 0 ; j < nOrigPts ; j++)
        {
            double pt[3];
            ifile >> pt[0];
            ifile >> pt[1];
            ifile >> pt[2];
            pts->SetPoint(i*nOrigPts+j, pt);
        }
    }

    int ncells = ds->GetNumberOfCells()*(nTimesteps-1);
    int ntris  = ds->GetNumberOfCells();
    ugrid->Allocate(ncells*(6+1));
    for (i = 0 ; i < ntris ; i++)
    {
        vtkCell *cell = ds->GetCell(i);
        if (cell->GetCellType() != VTK_TRIANGLE)
        {
            EXCEPTION1(InvalidVariableException, meshname);
        }
248
        vtkIdType wedge[6];
249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296
        wedge[0] = cell->GetPointId(0);
        wedge[1] = cell->GetPointId(1);
        wedge[2] = cell->GetPointId(2);
        for (j = 0 ; j < nTimesteps-1 ; j++)
        {
            wedge[3] = wedge[0] + nOrigPts;
            wedge[4] = wedge[1] + nOrigPts;
            wedge[5] = wedge[2] + nOrigPts;
            ugrid->InsertNextCell(VTK_WEDGE, 6, wedge);
            wedge[0] = wedge[3];
            wedge[1] = wedge[4];
            wedge[2] = wedge[5];
        }
    }
    ugrid->SetPoints(pts);
    pts->Delete();

    return ugrid;
}


// ****************************************************************************
//  Method: avtExtrudedVolFileFormat::GetVar
//
//  Purpose:
//      Gets a scalar variable associated with this file.  Although VTK has
//      support for many different types, the best bet is vtkFloatArray, since
//      that is supported everywhere through VisIt.
//
//  Arguments:
//      domain     The index of the domain.  If there are NDomains, this
//                 value is guaranteed to be between 0 and NDomains-1,
//                 regardless of block origin.
//      varname    The name of the variable requested.
//
//  Programmer: childs -- generated by xml2avt
//  Creation:   Fri May 18 17:52:04 PST 2007
//
// ****************************************************************************

vtkDataArray *
avtExtrudedVolFileFormat::GetVar(int domain, const char *varname)
{
    // We could probably do this more efficiently ... but, for now, read in
    // all the data one entry at a time.

    int  i, j, k;
    int  varId = -1;
297
    for (i = 0 ; i < (int)variables.size() ; i++)
298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367
    {
        if (variables[i] == varname)
            varId = i;
    }
    if (varId < 0)
        EXCEPTION1(InvalidVariableException, varname);

    vtkFloatArray *rv = vtkFloatArray::New();

    char filename[1024];
    for (i = 0 ; i < nTimesteps ; i++)
    {
        SNPRINTF(filename, 1024, "%s.%d.%d.exvol_var", stem.c_str(), i,domain);
        ifstream ifile(filename);
        if (ifile.fail())
        {
            EXCEPTION1(InvalidVariableException, varname);
        }
        int npts;
        ifile >> npts;
        if (i == 0)
            rv->SetNumberOfTuples(npts*nTimesteps);

        double val;
        for (k = 0 ; k < npts ; k++)
        {
            ifile >> val; // X
            ifile >> val; // Y
            ifile >> val; // Z
        }
        for (j = 0 ; j < varId ; j++)
            for (k = 0 ; k < npts ; k++)
                ifile >> val;
        for (k = 0 ; k < npts ; k++)
        {
            ifile >> val;
            rv->SetTuple1(i*npts+k, val);
        }
    }

    return rv;
}


// ****************************************************************************
//  Method: avtExtrudedVolFileFormat::GetVectorVar
//
//  Purpose:
//      Gets a vector variable associated with this file.  Although VTK has
//      support for many different types, the best bet is vtkFloatArray, since
//      that is supported everywhere through VisIt.
//
//  Arguments:
//      domain     The index of the domain.  If there are NDomains, this
//                 value is guaranteed to be between 0 and NDomains-1,
//                 regardless of block origin.
//      varname    The name of the variable requested.
//
//  Programmer: childs -- generated by xml2avt
//  Creation:   Fri May 18 17:52:04 PST 2007
//
// ****************************************************************************

vtkDataArray *
avtExtrudedVolFileFormat::GetVectorVar(int domain, const char *varname)
{
    EXCEPTION1(InvalidVariableException, varname);
}