File: FEDataStructures.cxx

package info (click to toggle)
paraview 5.11.0%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 497,236 kB
  • sloc: cpp: 3,171,290; ansic: 1,315,072; python: 134,290; xml: 103,324; sql: 65,887; sh: 5,286; javascript: 4,901; yacc: 4,383; java: 3,977; perl: 2,363; lex: 1,909; f90: 1,255; objc: 143; makefile: 119; tcl: 59; pascal: 50; fortran: 29
file content (159 lines) | stat: -rw-r--r-- 4,305 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
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
/*=========================================================================

  Program:   ParaView
  Module:    vtkDataObjectToConduit.h

  Copyright (c) Kitware, Inc.
  All rights reserved.
  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.

     This software is distributed WITHOUT ANY WARRANTY; without even
     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
     PURPOSE.  See the above copyright notice for more information.

=========================================================================*/

#include "FEDataStructures.h"

#include <iostream>
#include <iterator>
#include <mpi.h>

Grid::Grid() = default;

void Grid::Initialize(const unsigned int numPoints[3], const double spacing[3])
{
  if (numPoints[0] == 0 || numPoints[1] == 0 || numPoints[2] == 0)
  {
    std::cerr << "Must have a non-zero amount of points in each dimension.\n";
  }
  // in parallel, we do a simple partitioning in the x-direction.
  int mpiSize = 1;
  int mpiRank = 0;
  MPI_Comm_rank(MPI_COMM_WORLD, &mpiRank);
  MPI_Comm_size(MPI_COMM_WORLD, &mpiSize);

  unsigned int startXPoint = mpiRank * numPoints[0] / mpiSize;
  unsigned int endXPoint = (mpiRank + 1) * numPoints[0] / mpiSize;
  if (mpiSize != mpiRank + 1)
  {
    endXPoint++;
  }

  // create the points -- slowest in the x and fastest in the z directions
  double coord[3] = { 0, 0, 0 };
  for (unsigned int x = startXPoint; x < endXPoint; x++)
  {
    coord[0] = x * spacing[0];
    for (unsigned int y = 0; y < numPoints[1]; y++)
    {
      coord[1] = y * spacing[1];
      for (unsigned int z = 0; z < numPoints[2]; z++)
      {
        coord[2] = z * spacing[2];
        // add the coordinate to the end of the vector
        std::copy(coord, coord + 3, std::back_inserter(this->Points));
      }
    }
  }
  // create the hex cells
  unsigned int numXPoints = endXPoint - startXPoint;
  for (unsigned int i = 0; i < numXPoints - 1; i++)
  {
    for (unsigned int j = 0; j < numPoints[1] - 1; j++)
    {
      for (unsigned int k = 0; k < numPoints[2] - 1; k++)
      {
        unsigned int cellPoints[8] = { i * numPoints[1] * numPoints[2] + j * numPoints[2] + k,
          (i + 1) * numPoints[1] * numPoints[2] + j * numPoints[2] + k,
          (i + 1) * numPoints[1] * numPoints[2] + (j + 1) * numPoints[2] + k,
          i * numPoints[1] * numPoints[2] + (j + 1) * numPoints[2] + k,
          i * numPoints[1] * numPoints[2] + j * numPoints[2] + k + 1,
          (i + 1) * numPoints[1] * numPoints[2] + j * numPoints[2] + k + 1,
          (i + 1) * numPoints[1] * numPoints[2] + (j + 1) * numPoints[2] + k + 1,
          i * numPoints[1] * numPoints[2] + (j + 1) * numPoints[2] + k + 1 };
        std::copy(cellPoints, cellPoints + 8, std::back_inserter(this->Cells));
      }
    }
  }
}

size_t Grid::GetNumberOfPoints()
{
  return this->Points.size() / 3;
}

size_t Grid::GetNumberOfCells()
{
  return this->Cells.size() / 8;
}

double* Grid::GetPointsArray()
{
  if (this->Points.empty())
  {
    return nullptr;
  }
  return this->Points.data();
}

double* Grid::GetPoint(size_t pointId)
{
  if (pointId >= this->GetNumberOfPoints())
  {
    return nullptr;
  }
  return this->Points.data() + pointId * 3;
}

unsigned int* Grid::GetCellPoints(size_t cellId)
{
  if (cellId >= this->GetNumberOfCells())
  {
    return nullptr;
  }
  return this->Cells.data() + cellId * 8;
}

Attributes::Attributes()
{
  this->GridPtr = nullptr;
}

void Attributes::Initialize(Grid* grid)
{
  this->GridPtr = grid;
}

void Attributes::UpdateFields(double time)
{
  size_t numPoints = this->GridPtr->GetNumberOfPoints();
  this->Velocity.resize(numPoints * 3);
  for (size_t pt = 0; pt < numPoints; pt++)
  {
    double* coord = this->GridPtr->GetPoint(pt);
    this->Velocity[pt] = coord[1] * time;
  }
  std::fill(this->Velocity.begin() + numPoints, this->Velocity.end(), 0.);
  size_t numCells = this->GridPtr->GetNumberOfCells();
  this->Pressure.resize(numCells);
  std::fill(this->Pressure.begin(), this->Pressure.end(), 1.f);
}

double* Attributes::GetVelocityArray()
{
  if (this->Velocity.empty())
  {
    return nullptr;
  }
  return this->Velocity.data();
}

float* Attributes::GetPressureArray()
{
  if (this->Pressure.empty())
  {
    return nullptr;
  }
  return this->Pressure.data();
}