Chaste  Build::
AbstractRegularGridDiscreteContinuumSolver.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4  All rights reserved.
5 
6  University of Oxford means the Chancellor, Masters and Scholars of the
7  University of Oxford, having an administrative office at Wellington
8  Square, Oxford OX1 2JD, UK.
9 
10  This file is part of Chaste.
11 
12  Redistribution and use in source and binary forms, with or without
13  modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24  AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25  IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26  ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27  LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28  CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29  GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30  HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32  OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34  */
35 
36 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning for now (gcc4.3)
37 #include <vtkDoubleArray.h>
38 #include <vtkPointData.h>
39 #include <vtkPolyData.h>
40 #include <vtkProbeFilter.h>
41 #include <vtkImageData.h>
42 #include <vtkSmartPointer.h>
43 #include "RegularGridWriter.hpp"
44 
45 #include "AbstractRegularGridDiscreteContinuumSolver.hpp"
46 
47 template<unsigned DIM>
50  mpVtkSolution(),
51  mpRegularGrid()
52 {
53  this->mHasRegularGrid = true;
54 }
55 
56 template<unsigned DIM>
58 {
59 
60 }
61 
62 template<unsigned DIM>
63 boost::shared_ptr<RegularGrid<DIM> > AbstractRegularGridDiscreteContinuumSolver<DIM>::GetGrid()
64 {
65  if(!this->mpRegularGrid)
66  {
67  EXCEPTION("A regular grid has not been set.");
68  }
69  return this->mpRegularGrid;
70 }
71 
72 template<unsigned DIM>
73 std::vector<units::quantity<unit::concentration> > AbstractRegularGridDiscreteContinuumSolver<DIM>::GetConcentrations(const std::vector<DimensionalChastePoint<DIM> >& samplePoints)
74 {
75  if(!this->mpVtkSolution)
76  {
77  this->Setup();
78  }
79 
80  std::vector<units::quantity<unit::concentration> > sampled_solution(samplePoints.size(), 0.0*this->mReferenceConcentration);
81 
82  // Sample the field at these locations
83  vtkSmartPointer<vtkPolyData> p_polydata = vtkSmartPointer<vtkPolyData>::New();
84  vtkSmartPointer<vtkPoints> p_points = vtkSmartPointer<vtkPoints>::New();
85  p_points->SetNumberOfPoints(samplePoints.size());
86  for(unsigned idx=0; idx< samplePoints.size(); idx++)
87  {
88  c_vector<double, DIM> location = samplePoints[idx].GetLocation(this->mpRegularGrid->GetReferenceLengthScale());
89  if(DIM==3)
90  {
91  p_points->SetPoint(idx, location[0], location[1], location[2]);
92  }
93  else
94  {
95  p_points->SetPoint(idx, location[0], location[1], 0.0);
96  }
97  }
98  p_polydata->SetPoints(p_points);
99 
100  vtkSmartPointer<vtkProbeFilter> p_probe_filter = vtkSmartPointer<vtkProbeFilter>::New();
101  #if VTK_MAJOR_VERSION <= 5
102  p_probe_filter->SetInput(p_polydata);
103  p_probe_filter->SetSource(this->mpVtkSolution);
104  #else
105  p_probe_filter->SetInputData(p_polydata);
106  p_probe_filter->SetSourceData(this->mpVtkSolution);
107  #endif
108  p_probe_filter->Update();
109  vtkSmartPointer<vtkPointData> p_point_data = p_probe_filter->GetOutput()->GetPointData();
110 
111  unsigned num_points = p_point_data->GetArray(this->mLabel.c_str())->GetNumberOfTuples();
112  for(unsigned idx=0; idx<num_points; idx++)
113  {
114  sampled_solution[idx] = p_point_data->GetArray(this->mLabel.c_str())->GetTuple1(idx)*this->mReferenceConcentration;
115  }
116  return sampled_solution;
117 }
118 
119 template<unsigned DIM>
120 std::vector<units::quantity<unit::concentration> > AbstractRegularGridDiscreteContinuumSolver<DIM>::GetConcentrations(boost::shared_ptr<RegularGrid<DIM> > pGrid)
121 {
122  if(this->mpRegularGrid == pGrid)
123  {
124  return this->GetConcentrations();
125  }
126  else
127  {
128  return this->GetConcentrations(pGrid->GetLocations());
129  }
130 }
131 
132 template<unsigned DIM>
133 std::vector<units::quantity<unit::concentration> > AbstractRegularGridDiscreteContinuumSolver<DIM>::GetConcentrations(boost::shared_ptr<DiscreteContinuumMesh<DIM> > pMesh)
134 {
135  return this->GetConcentrations(pMesh->GetNodeLocationsAsPoints());
136 }
137 
138 template<unsigned DIM>
140 {
141  if(!this->mpVtkSolution)
142  {
143  this->Setup();
144  }
145 
146  std::vector<double> sampled_solution(rSamplePoints.size(), 0.0);
147 
148  // Sample the field at these locations
149  vtkSmartPointer<vtkPolyData> p_polydata = vtkSmartPointer<vtkPolyData>::New();
150  vtkSmartPointer<vtkPoints> p_points = vtkSmartPointer<vtkPoints>::New();
151  p_points->SetNumberOfPoints(rSamplePoints.size());
152  for(unsigned idx=0; idx< rSamplePoints.size(); idx++)
153  {
154  c_vector<double, DIM> location = rSamplePoints[idx].GetLocation(this->mpRegularGrid->GetReferenceLengthScale());
155  if(DIM==3)
156  {
157  p_points->SetPoint(idx, location[0], location[1], location[2]);
158  }
159  else
160  {
161  p_points->SetPoint(idx, location[0], location[1], 0.0);
162  }
163  }
164  p_polydata->SetPoints(p_points);
165 
166  vtkSmartPointer<vtkProbeFilter> p_probe_filter = vtkSmartPointer<vtkProbeFilter>::New();
167  #if VTK_MAJOR_VERSION <= 5
168  p_probe_filter->SetInput(p_polydata);
169  p_probe_filter->SetSource(this->mpVtkSolution);
170  #else
171  p_probe_filter->SetInputData(p_polydata);
172  p_probe_filter->SetSourceData(this->mpVtkSolution);
173  #endif
174  p_probe_filter->Update();
175  vtkSmartPointer<vtkPointData> p_point_data = p_probe_filter->GetOutput()->GetPointData();
176 
177  unsigned num_points = p_point_data->GetArray(this->mLabel.c_str())->GetNumberOfTuples();
178  for(unsigned idx=0; idx<num_points; idx++)
179  {
180  sampled_solution[idx] = p_point_data->GetArray(this->mLabel.c_str())->GetTuple1(idx);
181  }
182  return sampled_solution;
183 }
184 
185 template<unsigned DIM>
187 {
188  if(this->mpRegularGrid == pGrid)
189  {
190  return this->GetSolution();
191  }
192  else
193  {
194  return this->GetSolution(pGrid->GetLocations());
195  }
196 }
197 
198 template<unsigned DIM>
200 {
201  return this->GetSolution(pMesh->GetNodeLocationsAsPoints());
202 }
203 
204 template<unsigned DIM>
206 {
207  if(!this->mpVtkSolution)
208  {
209  this->Setup();
210  }
211  return this->mpVtkSolution;
212 }
213 
214 template<unsigned DIM>
216 {
217  this->mpRegularGrid = pGrid;
218 }
219 
220 template<unsigned DIM>
222 {
223  if(!this->mpRegularGrid)
224  {
225  EXCEPTION("Regular grid DiscreteContinuum solvers need a grid before Setup can be called.");
226  }
227 
228  // Set up the VTK solution
229  this->mpVtkSolution = vtkSmartPointer<vtkImageData>::New();
230 
231  if(DIM==3)
232  {
233  this->mpVtkSolution->SetDimensions(this->mpRegularGrid->GetExtents()[0], this->mpRegularGrid->GetExtents()[1], this->mpRegularGrid->GetExtents()[2]);
234  }
235  else
236  {
237  this->mpVtkSolution->SetDimensions(this->mpRegularGrid->GetExtents()[0], this->mpRegularGrid->GetExtents()[1], 1);
238  }
239 
240  double spacing = this->mpRegularGrid->GetSpacing()/this->mpRegularGrid->GetReferenceLengthScale();
241  this->mpVtkSolution->SetSpacing(spacing, spacing, spacing);
242 
243  c_vector<double,DIM> origin = this->mpRegularGrid->GetOrigin().GetLocation(this->mpRegularGrid->GetReferenceLengthScale());
244  if(DIM==3)
245  {
246  this->mpVtkSolution->SetOrigin(origin[0], origin[1], origin[2]);
247  }
248  else
249  {
250  this->mpVtkSolution->SetOrigin(origin[0], origin[1], 0.0);
251  }
252 
253  this->mSolution = std::vector<double>(0.0, this->mpRegularGrid->GetNumberOfPoints());
254 }
255 
256 template<unsigned DIM>
258 {
259  if(!this->mpVtkSolution)
260  {
261  this->Setup();
262  }
263 
264  vtkSmartPointer<vtkDoubleArray> pPointData = vtkSmartPointer<vtkDoubleArray>::New();
265  pPointData->SetNumberOfComponents(1);
266  pPointData->SetNumberOfTuples(data.size());
267  pPointData->SetName(this->GetLabel().c_str());
268  for (unsigned i = 0; i < data.size(); i++)
269  {
270  pPointData->SetValue(i, data[i]);
271  }
272  this->mpVtkSolution->GetPointData()->AddArray(pPointData);
273 
274  // Note, if the data vector being passed in is mPointSolution, then it will be overwritten with zeros.
275  this->mSolution = std::vector<double>(data.size(), 0.0);
276  for (unsigned i = 0; i < data.size(); i++)
277  {
278  this->mSolution[i] = data[i];
279  }
280 
281  this->mConcentrations = std::vector<units::quantity<unit::concentration> >(data.size(), 0.0*this->mReferenceConcentration);
282  for (unsigned i = 0; i < data.size(); i++)
283  {
284  this->mConcentrations[i] = data[i]*this->mReferenceConcentration;
285  }
286 }
287 
288 template<unsigned DIM>
289 void AbstractRegularGridDiscreteContinuumSolver<DIM>::UpdateSolution(std::vector<units::quantity<unit::concentration> >& data)
290 {
291  if(!this->mpVtkSolution)
292  {
293  this->Setup();
294  }
295 
296  vtkSmartPointer<vtkDoubleArray> pPointData = vtkSmartPointer<vtkDoubleArray>::New();
297  pPointData->SetNumberOfComponents(1);
298  pPointData->SetNumberOfTuples(data.size());
299  pPointData->SetName((this->GetLabel()).c_str());
300  for (unsigned i = 0; i < data.size(); i++)
301  {
302  pPointData->SetValue(i, data[i]/this->mReferenceConcentration);
303  }
304  this->mpVtkSolution->GetPointData()->AddArray(pPointData);
305 
306  // Note, if the data vector being passed in is mPointSolution, then it will be overwritten with zeros.
307  this->mConcentrations = std::vector<units::quantity<unit::concentration> >(data.size(), 0.0*this->mReferenceConcentration);
308  for (unsigned i = 0; i < data.size(); i++)
309  {
310  this->mConcentrations[i] = data[i];
311  }
312 }
313 
314 template<unsigned DIM>
316 {
317  if(!this->mpVtkSolution)
318  {
319  this->Setup();
320  }
321 
322  if(!this->CellPopulationIsSet())
323  {
324  EXCEPTION("The DiscreteContinuum solver needs a cell population for this operation.");
325  }
326 
327  this->mpRegularGrid->SetCellPopulation(*(this->mpCellPopulation), this->mCellPopulationReferenceLength);
328  std::vector<std::vector<CellPtr> > point_cell_map = this->mpRegularGrid->GetPointCellMap();
329  for(unsigned idx=0; idx<point_cell_map.size(); idx++)
330  {
331  for(unsigned jdx=0; jdx<point_cell_map[idx].size(); jdx++)
332  {
333  point_cell_map[idx][jdx]->GetCellData()->SetItem(this->mLabel, this->mConcentrations[idx]/this->mCellPopulationReferenceConcentration);
334  }
335  }
336 }
337 
338 template<unsigned DIM>
340 {
341 
342 }
343 
344 template<unsigned DIM>
346 {
347  if(!this->mpVtkSolution)
348  {
349  this->Setup();
350  }
351 
352  if(!this->mpOutputFileHandler)
353  {
354  EXCEPTION("An output file handler has not been set for the DiscreteContinuum solver.");
355  }
356 
357  RegularGridWriter writer;
358  if(!this->mFilename.empty())
359  {
360  writer.SetFilename((this->mpOutputFileHandler->GetOutputDirectoryFullPath() + "/" + this->mFilename+".vti"));
361  }
362  else
363  {
364  writer.SetFilename((this->mpOutputFileHandler->GetOutputDirectoryFullPath() + "/solution.vti"));
365  }
366  writer.SetImage(this->mpVtkSolution);
367  writer.Write();
368 }
369 
370 // Explicit instantiation
virtual std::vector< units::quantity< unit::concentration > > GetConcentrations()
void SetFilename(const std::string &rFilename)
boost::shared_ptr< OutputFileHandler > mpOutputFileHandler
AbstractCellPopulation< DIM > * mpCellPopulation
std::vector< units::quantity< unit::concentration > > mConcentrations
units::quantity< unit::concentration > mReferenceConcentration
void SetImage(vtkSmartPointer< vtkImageData > pImage)
units::quantity< unit::concentration > mCellPopulationReferenceConcentration
void SetGrid(boost::shared_ptr< RegularGrid< DIM > > pRegularGrid)
units::quantity< unit::length > mCellPopulationReferenceLength