Scippy

UG

Ubiquity Generator framework

scipParaSolutionMpi.cpp
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the program and software framework */
4/* UG --- Ubquity Generator Framework */
5/* */
6/* Copyright Written by Yuji Shinano <shinano@zib.de>, */
7/* Copyright (C) 2021-2024 by Zuse Institute Berlin, */
8/* licensed under LGPL version 3 or later. */
9/* Commercial licenses are available through <licenses@zib.de> */
10/* */
11/* This code is free software; you can redistribute it and/or */
12/* modify it under the terms of the GNU Lesser General Public License */
13/* as published by the Free Software Foundation; either version 3 */
14/* of the License, or (at your option) any later version. */
15/* */
16/* This program is distributed in the hope that it will be useful, */
17/* but WITHOUT ANY WARRANTY; without even the implied warranty of */
18/* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
19/* GNU Lesser General Public License for more details. */
20/* */
21/* You should have received a copy of the GNU Lesser General Public License */
22/* along with this program. If not, see <http://www.gnu.org/licenses/>. */
23/* */
24/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
25
26/**@file scipParaSolutionMpi.cpp
27 * @brief ScipParaSolution extension for MPI communication.
28 * @author Yuji Shinano
29 *
30 *
31 *
32 */
33
34/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
35
36
37#include <mpi.h>
38#include "scipParaCommMpi.h"
39#include "scipParaSolutionMpi.h"
40
41using namespace UG;
42using namespace ParaSCIP;
43
44/** create clone of this object */
46ScipParaSolutionMpi::clone(ParaComm *comm)
47{
49}
50
51/** create ScipDiffSubproblemPreDatatype */
52MPI_Datatype
54 )
55{
56 const int nBlocks = 2;
57 MPI_Datatype preDatatype;
58
59 MPI_Aint startAddress = 0;
60 MPI_Aint address = 0;
61
62 int blockLengths[nBlocks];
63 MPI_Aint displacements[nBlocks];
64 MPI_Datatype types[nBlocks];
65
66 for( int i = 0; i < nBlocks; i++ ){
67 blockLengths[i] = 1;
68 }
69
71 MPI_Get_address( &objectiveFunctionValue, &startAddress )
72 );
73 displacements[0] = 0;
74 types[0] = MPI_DOUBLE;
75
77 MPI_Get_address( &nVars, &address )
78 );
79 displacements[1] = address - startAddress;
80 types[1] = MPI_INT;
81
83 MPI_Type_create_struct(nBlocks, blockLengths, displacements, types, &preDatatype)
84 );
85
86 return preDatatype;
87}
88
89/** create ScipDiffSubproblemDatatype */
90MPI_Datatype
92 bool memAllocNecessary
93 )
94{
95 const int nBlocks = 2;
96
97 MPI_Datatype datatype = MPI_DATATYPE_NULL;
98
99 MPI_Aint startAddress = 0;
100 MPI_Aint address = 0;
101
102 int blockLengths[nBlocks];
103 MPI_Aint displacements[nBlocks];
104 MPI_Datatype types[nBlocks];
105
106 if( nVars )
107 {
108 if( memAllocNecessary )
109 {
110 indicesAmongSolvers = new int[nVars];
111 values = new SCIP_Real[nVars];
112 }
113
114 MPI_CALL(
115 MPI_Get_address( indicesAmongSolvers, &startAddress )
116 );
117 displacements[0] = 0;
118 blockLengths[0] = nVars;
119 types[0] = MPI_INT;
120
121 MPI_CALL(
122 MPI_Get_address( values, &address )
123 );
124 displacements[1] = address - startAddress;
125 blockLengths[1] = nVars;
126 types[1] = MPI_DOUBLE;
127
128 MPI_CALL(
129 MPI_Type_create_struct(nBlocks, blockLengths, displacements, types, &datatype)
130 );
131 }
132 return datatype;
133}
134
135/** send solution data to the rank */
136void
138{
139 DEF_PARA_COMM( commMpi, comm);
140 MPI_Datatype preDatatype;
141 preDatatype = createPreDatatype();
142 MPI_CALL(
143 MPI_Type_commit( &preDatatype )
144 );
146 commMpi->ubcast(&objectiveFunctionValue, 1, preDatatype, root)
147 );
148 MPI_CALL(
149 MPI_Type_free( &preDatatype )
150 );
151
152 if( nVars ){
153 MPI_Datatype datatype;
154 if( comm->getRank() == root )
155 {
156 datatype = createDatatype(false);
157 }
158 else
159 {
160 datatype = createDatatype(true);
161 }
162 MPI_CALL(
163 MPI_Type_commit( &datatype )
164 );
166 commMpi->ubcast(indicesAmongSolvers, 1, datatype, root)
167 );
168 MPI_CALL(
169 MPI_Type_free( &datatype )
170 );
171 }
172}
173
174
175/** send solution data to the rank */
176void
178{
179 DEF_PARA_COMM( commMpi, comm);
180 MPI_Datatype preDatatype;
181 preDatatype = createPreDatatype();
182 MPI_CALL(
183 MPI_Type_commit( &preDatatype )
184 );
186 commMpi->usend(&objectiveFunctionValue, 1, preDatatype, destination, TagSolution)
187 );
188 MPI_CALL(
189 MPI_Type_free( &preDatatype )
190 );
191
192 if( nVars ){
193 MPI_Datatype datatype;
194 datatype = createDatatype(false);
195 MPI_CALL(
196 MPI_Type_commit( &datatype )
197 );
199 commMpi->usend(indicesAmongSolvers, 1, datatype, destination, TagSolution1)
200 );
201 MPI_CALL(
202 MPI_Type_free( &datatype )
203 );
204 }
205}
206
207/** receive solution data from the source rank */
208void
210{
211 DEF_PARA_COMM( commMpi, comm);
212 MPI_Datatype preDatatype;
213 preDatatype = createPreDatatype();
214 MPI_CALL(
215 MPI_Type_commit( &preDatatype )
216 );
218 commMpi->ureceive(&objectiveFunctionValue, 1, preDatatype, source, TagSolution)
219 );
220 MPI_CALL(
221 MPI_Type_free( &preDatatype )
222 );
223
224 if( nVars ){
225 MPI_Datatype datatype;
226 datatype = createDatatype(true);
227 MPI_CALL(
228 MPI_Type_commit( &datatype )
229 );
230 MPI_CALL(
231 MPI_Type_commit( &datatype )
232 );
234 commMpi->ureceive(indicesAmongSolvers, 1, datatype, source, TagSolution1)
235 );
236 MPI_CALL(
237 MPI_Type_free( &datatype )
238 );
239 }
240}
MPI_Datatype createDatatype(bool memAllocNecessary)
void send(UG::ParaComm *comm, int destination)
void bcast(UG::ParaComm *comm, int root)
void receive(UG::ParaComm *comm, int source)
int getRank()
get rank of caller's thread
Base class of communicator object.
Definition: paraComm.h:102
static ScipParaCommTh * comm
Definition: fscip.cpp:73
static const int TagSolution1
static const int TagSolution
Definition: paraTagDef.h:51
#define DEF_PARA_COMM(para_comm, comm)
#define MPI_CALL(mpicall)
Definition: paraCommMpi.h:68
#define PARA_COMM_CALL(paracommcall)
Definition: paraComm.h:47
SCIP ParaComm extension for MPI communication.
ScipParaSolution extension for MPI communication.