Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
ParallelColumnDataWriter.cpp
1/*
2
3Copyright (c) 2005-2024, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, 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
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#include "ParallelColumnDataWriter.hpp"
37#include "Exception.hpp"
38#include "DistributedVectorFactory.hpp"
39
41 const std::string& rBaseName,
42 bool cleanDirectory)
43 : ColumnDataWriter(rDirectory, rBaseName, cleanDirectory),
44 mConcentrated(nullptr)
45{
46 int num_procs;
47 MPI_Comm_size(PETSC_COMM_WORLD, &num_procs);
48 if (num_procs==1)
49 {
50 mIsParallel = false;
51 }
52 else
53 {
54 mIsParallel = true;
55 }
56}
57
58void ParallelColumnDataWriter::PutVector(int variableID, Vec petscVector)
59{
60 int size;
61 VecGetSize(petscVector,&size);
62
63 if (size != mFixedDimensionSize)
64 {
65 EXCEPTION("Size of vector does not match FixedDimensionSize.");
66 }
67
68 // Construct the appropriate "scatter" object to concentrate the vector on the master
69 if (mConcentrated==nullptr)
70 {
71 VecScatterCreateToZero(petscVector, &mToMaster, &mConcentrated);
72 }
73
74// int size2;
75// VecGetSize(mConcentrated, &size2);
76// std::cout << "Vector size=" << size << "," << size2 << std::endl << std::flush;
77
78//PETSc-3.x.x or PETSc-2.3.3
79#if ((PETSC_VERSION_MAJOR == 3) || (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 3)) //2.3.3 or 3.x.x
80 VecScatterBegin(mToMaster, petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD);
81 VecScatterEnd(mToMaster, petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD);
82#else
83 VecScatterBegin(petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD, mToMaster);
84 VecScatterEnd(petscVector, mConcentrated, INSERT_VALUES, SCATTER_FORWARD, mToMaster);
85#endif
86
87// std::cout << "Done scatter" << std::endl << std::flush;
88
90 {
91 double *concentrated_vector;
92 VecGetArray(mConcentrated, &concentrated_vector);
93 for (int i=0; i<size; i++)
94 {
95 ColumnDataWriter::PutVariable(variableID, concentrated_vector[i], i);
96 }
97 VecRestoreArray(mConcentrated, &concentrated_vector);
98 }
99}
100
102{
103 // Put the stripe into its own 'unstriped' vector
104 DistributedVectorFactory* p_factory = rStripe.GetFactory();
105 Vec unstriped_petsc = p_factory->CreateVec();
106 DistributedVector unstriped = p_factory->CreateDistributedVector(unstriped_petsc);
107 for (DistributedVector::Iterator index = unstriped.Begin();
108 index!= unstriped.End();
109 ++index)
110 {
111 unstriped[index] = rStripe[index];
112 }
113
114 // Put the unstriped vector
115 ParallelColumnDataWriter::PutVector(variableId, unstriped_petsc);
116 PetscTools::Destroy(unstriped_petsc);
117}
118
120{
122 {
124 }
125 else
126 {
127 mIsInDefineMode = false;
128 }
129}
130
139void ParallelColumnDataWriter::PutVariable(int variableID, double variableValue, long dimensionPosition)
140{
142 {
143 // Master process is allowed to write
144 ColumnDataWriter::PutVariable(variableID, variableValue, dimensionPosition);
145 }
146}
147
149{
150 if (mConcentrated != nullptr)
151 {
152 VecScatterDestroy(PETSC_DESTROY_PARAM(mToMaster));
154 }
155 Close();
156}
157
159{
160 // Paranoia
161 PetscTools::Barrier("ParallelColumnDataWriter::AdvanceAlongUnlimitedDimension");
162
164 {
166 }
167}
168
170{
171 // Paranoia
172 PetscTools::Barrier("ParallelColumnDataWriter::Close");
173
175 {
177 }
178}
#define EXCEPTION(message)
#define PETSC_DESTROY_PARAM(x)
virtual void EndDefineMode()
void DoAdvanceAlongUnlimitedDimension()
virtual void PutVariable(int variableID, double variableValue, long dimensionPosition=-1)
DistributedVector CreateDistributedVector(Vec vec, bool readOnly=false)
DistributedVectorFactory * GetFactory()
ParallelColumnDataWriter(const std::string &rDirectory, const std::string &rBaseName, bool cleanDirectory=true)
void PutVectorStripe(int variableId, DistributedVector::Stripe &rStripe)
void PutVariable(int variableID, double variableValue, long dimensionPosition=-1)
void PutVector(int variableID, Vec petscVector)
static void Destroy(Vec &rVec)
static bool AmMaster()
static void Barrier(const std::string callerId="")