GeneratedFieldPackInfoTestGPU.cpp 7.48 KB
Newer Older
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
160
161
162
163
164
165
166
167
168
169
170
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
205
206
207
208
209
210
211
212
213
//======================================================================================================================
//
//  This file is part of waLBerla. waLBerla is free software: you can 
//  redistribute it and/or modify it under the terms of the GNU General Public
//  License as published by the Free Software Foundation, either version 3 of 
//  the License, or (at your option) any later version.
//  
//  waLBerla is distributed in the hope that it will be useful, but WITHOUT 
//  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 
//  FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License 
//  for more details.
//  
//  You should have received a copy of the GNU General Public License along
//  with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
//! \file GeneratedFieldPackInfoTestGPU.cpp
//! \ingroup field
//! \author Helen Schottenhamml <helen.schottenhamml@fau.de>
//! \brief Tests if a GPU Field is correctly communicated using generated pack info
//
//======================================================================================================================

#include "field/AddToStorage.h"
#include "field/GhostLayerField.h"

#include "blockforest/Initialization.h"

#include "core/debug/TestSubsystem.h"
#include "core/Environment.h"

#include "cuda/FieldCopy.h"
#include "cuda/communication/UniformGPUScheme.h"

#include "stencil/D3Q27.h"

// include generated files
#include "ScalarFieldCommunicationGPU.h"
#include "ScalarFieldPullReductionGPU.h"


namespace walberla {

using Stencil_T = stencil::D3Q27;

cuda::GPUField<int> * createGPUField( IBlock* const block, StructuredBlockStorage* const storage ) {

   return new cuda::GPUField<int> (
      storage->getNumberOfXCells( *block ), // number of cells in x direction
      storage->getNumberOfYCells( *block ), // number of cells in y direction
      storage->getNumberOfZCells( *block ), // number of cells in z direction
      1,                                    // fSize
      1,                                    // number of ghost layers
      field::fzyx );

}

cuda::GPUField<int> * createSmallGPUField( IBlock * const , StructuredBlockStorage * const ) {
   return new cuda::GPUField<int> (2, 2, 2, 1, 1, field::fzyx );
}


void testScalarField( std::shared_ptr<blockforest::StructuredBlockForest> & sbf, BlockDataID gpuFieldId ) {

   cuda::communication::UniformGPUScheme< Stencil_T > us{ sbf };
   us.addPackInfo(std::make_shared< pystencils::ScalarFieldCommunicationGPU >(gpuFieldId));

   for( auto & block : *sbf ) {

      auto & gpuField = *(block.getData< cuda::GPUField< int > >(gpuFieldId));

      field::GhostLayerField< int, 1 > cpuField(gpuField.xSize(), gpuField.ySize(), gpuField.zSize(), 1, 0,
                                                field::fzyx);
      cpuField.setWithGhostLayer(0);

      WALBERLA_CHECK_EQUAL(cpuField.xSize(), 2)
      WALBERLA_CHECK_EQUAL(cpuField.ySize(), 2)
      WALBERLA_CHECK_EQUAL(cpuField.zSize(), 2)

      // initialize the bottom boundary
      cpuField(0, 0, 0) = 1;
      cpuField(0, 1, 0) = 2;
      cpuField(1, 0, 0) = 3;
      cpuField(1, 1, 0) = 4;

      cuda::fieldCpy(gpuField, cpuField);

      // communicate
      us.communicate();

      cuda::fieldCpy(cpuField, gpuField);

      WALBERLA_CHECK_EQUAL(cpuField(0, 0, +2), 1)
      WALBERLA_CHECK_EQUAL(cpuField(0, 1, +2), 2)
      WALBERLA_CHECK_EQUAL(cpuField(1, 0, +2), 3)
      WALBERLA_CHECK_EQUAL(cpuField(1, 1, +2), 4)

   }
}

void testScalarFieldPullReduction( std::shared_ptr<blockforest::StructuredBlockForest> & sbf, BlockDataID gpuFieldId ) {

   cuda::communication::UniformGPUScheme< Stencil_T > us1{ sbf };
   us1.addPackInfo(std::make_shared< pystencils::ScalarFieldPullReductionGPU >(gpuFieldId));

   cuda::communication::UniformGPUScheme< Stencil_T > us2{ sbf };
   us2.addPackInfo(std::make_shared< pystencils::ScalarFieldCommunicationGPU >(gpuFieldId));

   for( auto & block : *sbf ) {

      auto& gpuField = *(block.getData< cuda::GPUField< int > >(gpuFieldId));

      field::GhostLayerField< int, 1 > cpuField(gpuField.xSize(), gpuField.ySize(), gpuField.zSize(), 1, 0,
                                                field::fzyx);
      cpuField.setWithGhostLayer(0);

      WALBERLA_CHECK_EQUAL(cpuField.xSize(), 2)
      WALBERLA_CHECK_EQUAL(cpuField.ySize(), 2)
      WALBERLA_CHECK_EQUAL(cpuField.zSize(), 2)

      // initialize the bottom ghost layer cells
      cpuField(0, 0, -1) = 1;
      cpuField(0, 1, -1) = 2;
      cpuField(1, 0, -1) = 3;
      cpuField(1, 1, -1) = 4;

      // initialize the top interior cells
      cpuField(0, 0, 1) = 1;
      cpuField(0, 1, 1) = 1;
      cpuField(1, 0, 1) = 1;
      cpuField(1, 1, 1) = 1;

      cuda::fieldCpy(gpuField, cpuField);

      // communicate pull += reduction
      us1.communicate();

      cuda::fieldCpy(cpuField, gpuField);

      // check values in top ghost layer
      WALBERLA_CHECK_EQUAL(cpuField(0, 0, 2), 0)
      WALBERLA_CHECK_EQUAL(cpuField(0, 1, 2), 0)
      WALBERLA_CHECK_EQUAL(cpuField(1, 0, 2), 0)
      WALBERLA_CHECK_EQUAL(cpuField(1, 1, 2), 0)

      // check values in top interior cells
      WALBERLA_CHECK_EQUAL(cpuField(0, 0, 1), 2)
      WALBERLA_CHECK_EQUAL(cpuField(0, 1, 1), 3)
      WALBERLA_CHECK_EQUAL(cpuField(1, 0, 1), 4)
      WALBERLA_CHECK_EQUAL(cpuField(1, 1, 1), 5)

      // communicate to sync ghost layers
      us2.communicate();

      cuda::fieldCpy(cpuField, gpuField);

      // check values in bottom ghost layer
      WALBERLA_CHECK_EQUAL(cpuField(0, 0, -1), 2)
      WALBERLA_CHECK_EQUAL(cpuField(0, 1, -1), 3)
      WALBERLA_CHECK_EQUAL(cpuField(1, 0, -1), 4)
      WALBERLA_CHECK_EQUAL(cpuField(1, 1, -1), 5)

      // check values in top interior cells
      WALBERLA_CHECK_EQUAL(cpuField(0, 0, 1), 2)
      WALBERLA_CHECK_EQUAL(cpuField(0, 1, 1), 3)
      WALBERLA_CHECK_EQUAL(cpuField(1, 0, 1), 4)
      WALBERLA_CHECK_EQUAL(cpuField(1, 1, 1), 5)

   }
}

int main(int argc, char **argv) {

   using blockforest::createUniformBlockGrid;

   debug::enterTestMode();
   Environment walberlaEnv(argc,argv);

   // Create a BlockForest with 2x2x2 cells per block
   uint_t processes = uint_c( MPIManager::instance()->numProcesses() );
   auto blocks = createUniformBlockGrid(processes,1 ,1, //blocks
                                        2,2,2, //cells
                                        1, //dx
                                        true, //one block per process
                                        true,true,true); //periodicity

   // Create a Field with the same number of cells as the block
   BlockDataID scalarGPUFieldId = blocks->addStructuredBlockData<cuda::GPUField<int> > ( &createGPUField, "ScalarGPUField" );

   testScalarField( blocks, scalarGPUFieldId );

   // Create a BlockForest with 8x8x8 cells per block
   blocks = createUniformBlockGrid(processes,1 ,1, //blocks
                                   8,8,8,          //cells
                                   1,              //dx
                                   true,          //one block per process
                                   true,true,true);//periodicity

   // Create a Field with one quarter as many cells per dimension, i.e. a field with the same size as the one above
   scalarGPUFieldId = blocks->addStructuredBlockData<cuda::GPUField<int> > ( &createSmallGPUField, "ScalarGPUField" );

   testScalarField( blocks, scalarGPUFieldId );

   testScalarFieldPullReduction( blocks, scalarGPUFieldId );

   return 0;

}

} // namespace walberla

int main( int argc, char* argv[] ) {
   return walberla::main( argc, argv );
}