00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035 #include "ActorPetscSOE.h"
00036 #include "PetscSolver.h"
00037 #include "PetscSOE.h"
00038 #include <Matrix.h>
00039 #include <Graph.h>
00040 #include <Vertex.h>
00041 #include <VertexIter.h>
00042 #include <f2c.h>
00043 #include <math.h>
00044 #include <Channel.h>
00045 #include <FEM_ObjectBroker.h>
00046
00047 ActorPetscSOE::ActorPetscSOE(PetscSolver &theSOESolver, int blockSize)
00048 :myRank(0), theSOE(0), theSolver(&theSOESolver), recvBuffer(0)
00049 {
00050
00051 MPI_Comm_rank(PETSC_COMM_WORLD, &myRank);
00052 MPI_Comm_size(PETSC_COMM_WORLD, &numProcessors);
00053
00054 if (myRank == 0) {
00055 opserr << " ActorPetscSOE::ActorPetscSOE - must be rank 0\n";
00056 }
00057 recvBuffer = (void *)(&recvData[0]);
00058 MPI_Barrier(PETSC_COMM_WORLD);
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130 theSOE = new PetscSOE(*theSolver, blockSize);
00131 }
00132
00133
00134
00135 ActorPetscSOE::~ActorPetscSOE()
00136 {
00137
00138
00139 }
00140
00141
00142 int
00143 ActorPetscSOE::run(void)
00144 {
00145 int flag = 1;
00146 int tag, low, high, ierr, n, numRows;
00147 double *theData;
00148 int dnz = 0;
00149 int onz = 0;
00150 int *dnnz = PETSC_NULL;
00151 int *onnz = PETSC_NULL;
00152 MPI_Status status;
00153 void *buffer = 0;
00154
00155 while (flag != 0) {
00156 MPI_Bcast(recvBuffer, 3, MPI_INT, 0, PETSC_COMM_WORLD);
00157 flag = recvData[0];
00158 switch(flag) {
00159
00160 case 0:
00161 MPI_Barrier(PETSC_COMM_WORLD);
00162 break;
00163
00164 case 1:
00165 if (theSOE == 0)
00166 return -1;
00167 theSOE->isFactored = recvData[1];
00168 theSOE->solve();
00169 break;
00170
00171 case 2:
00172 if (theSOE == 0)
00173 return -1;
00174 tag = 100;
00175 MPI_Recv(recvBuffer, 3, MPI_INT, 0, tag, PETSC_COMM_WORLD, &status);
00176
00177 numRows = recvData[1];
00178 n = recvData[2];
00179 onnz = new int[numRows];
00180 dnnz = new int[numRows];
00181 buffer = (void *)dnnz;
00182 tag = 101;
00183 MPI_Recv(buffer, numRows, MPI_INT, 0, tag, PETSC_COMM_WORLD, &status);
00184
00185 buffer = (void *)onnz;
00186 tag = 102;
00187 MPI_Recv(buffer, numRows, MPI_INT, 0, tag, PETSC_COMM_WORLD, &status);
00188
00189 MPI_Bcast(recvBuffer, 3, MPI_INT, 0, PETSC_COMM_WORLD);
00190 theSOE->setSizeParallel(numRows, n, dnz, dnnz, onz, onnz);
00191 break;
00192
00193 case 3:
00194 if (theSOE == 0)
00195 return -1;
00196 theSOE->zeroA();
00197 break;
00198
00199 case 4:
00200 if (theSOE == 0)
00201 return -1;
00202 theSOE->zeroB();
00203 break;
00204
00205 case 5:
00206 if (theSOE == 0)
00207 return -1;
00208 tag = 99;
00209 ierr = VecGetOwnershipRange(theSOE->x, &low, &high); CHKERRA(ierr);
00210 recvData[0] = low; recvData[1] = high;
00211 MPI_Send(recvBuffer, 2, MPI_INT, 0, tag, PETSC_COMM_WORLD);
00212 ierr = VecGetArray(theSOE->x, &theData); CHKERRA(ierr);
00213 MPI_Send(theData, high-low, MPI_DOUBLE, 0, tag, PETSC_COMM_WORLD);
00214 ierr = VecRestoreArray(theSOE->x, &theData); CHKERRA(ierr);
00215 break;
00216
00217 case 6:
00218 if (theSOE == 0)
00219 return -1;
00220
00221
00222 default:
00223 opserr << "ActorPetscSOE::invalid action " << flag << " received\n";
00224 }
00225 }
00226 return 0;
00227 }
00228
00229
00230
00231