00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #ifdef VMDMPI
00018 #include <mpi.h>
00019 #include <stdio.h>
00020 #include <string.h>
00021 #include <Inform.h>
00022 #include <utilities.h>
00023 #include <WKFThreads.h>
00024 #include "VMDMPI.h"
00025
00026 typedef struct {
00027 int numcpus;
00028 int numgpus;
00029 float cpuspeed;
00030 float nodespeed;
00031 char machname[512];
00032 long corefree;
00033 long corepcnt;
00034 } nodeinfo;
00035
00036
00037 int vmd_mpi_init(int *argc, char ***argv) {
00038 int numnodes=0, mynode=0;
00039
00040 #if defined(VMDTHREADS)
00041 int provided=0;
00042 MPI_Init_thread(argc, argv, MPI_THREAD_SERIALIZED, &provided);
00043 if (provided != MPI_THREAD_SERIALIZED) {
00044 msgWarn << "MPI not providing thread-serial access." << sendmsg;
00045 }
00046 #else
00047 MPI_Init(argc, argv);
00048 #endif
00049 MPI_Comm_rank(MPI_COMM_WORLD, &mynode);
00050 MPI_Comm_size(MPI_COMM_WORLD, &numnodes);
00051
00052
00053 if (mynode != 0) {
00054 msgInfo.mute();
00055 msgWarn.mute();
00056 msgErr.mute();
00057 }
00058
00059 return 0;
00060 }
00061
00062 int vmd_mpi_barrier() {
00063 MPI_Barrier(MPI_COMM_WORLD);
00064 return 0;
00065 }
00066
00067 int vmd_mpi_fini() {
00068 vmd_mpi_barrier();
00069 msgInfo << "All nodes have reached the MPI shutdown barrier." << sendmsg;
00070
00071 MPI_Finalize();
00072
00073 return 0;
00074 }
00075
00076 int vmd_mpi_nodescan(int *noderank, int *nodecount,
00077 char *nodename, int maxnodenamelen,
00078 int localgpucount) {
00079 int numnodes=0, mynode=0;
00080 int namelen;
00081 char namebuf[MPI_MAX_PROCESSOR_NAME];
00082
00083 MPI_Comm_rank(MPI_COMM_WORLD, &mynode);
00084 MPI_Comm_size(MPI_COMM_WORLD, &numnodes);
00085 msgInfo << "Initializing parallel VMD instances via MPI..." << sendmsg;
00086
00087 nodeinfo *nodes;
00088 nodes = (nodeinfo *) malloc(numnodes * sizeof(nodeinfo));
00089 nodes[mynode].numcpus = wkf_thread_numprocessors();
00090 nodes[mynode].numgpus = localgpucount;
00091 nodes[mynode].cpuspeed = 1.0;
00092 nodes[mynode].nodespeed = nodes[mynode].numcpus * nodes[mynode].cpuspeed;
00093 nodes[mynode].corefree = vmd_get_avail_physmem_mb();
00094 nodes[mynode].corepcnt = vmd_get_avail_physmem_percent();
00095
00096 MPI_Get_processor_name((char *) &namebuf, &namelen);
00097
00098
00099 strncpy((char *) &nodes[mynode].machname, namebuf,
00100 (((namelen + 1) < 511) ? (namelen+1) : 511));
00101
00102
00103 strncpy(nodename, namebuf,
00104 (((namelen + 1) < (maxnodenamelen-1)) ? (namelen+1) : (maxnodenamelen-1)));
00105
00106 #if defined(USE_MPI_IN_PLACE)
00107
00108 MPI_Allgather(MPI_IN_PLACE, sizeof(nodeinfo), MPI_BYTE,
00109 &nodes[ 0], sizeof(nodeinfo), MPI_BYTE,
00110 MPI_COMM_WORLD);
00111 #else
00112
00113 MPI_Allgather(&nodes[mynode], sizeof(nodeinfo), MPI_BYTE,
00114 &nodes[ 0], sizeof(nodeinfo), MPI_BYTE,
00115 MPI_COMM_WORLD);
00116 #endif
00117
00118 if (mynode == 0) {
00119 char msgtxt[1024];
00120 float totalspeed = 0.0;
00121 int totalcpus = 0;
00122 int totalgpus = 0;
00123 int i;
00124
00125 for (i=0; i<numnodes; i++) {
00126 totalcpus += nodes[i].numcpus;
00127 totalgpus += nodes[i].numgpus;
00128 totalspeed += nodes[i].nodespeed;
00129 }
00130
00131 sprintf(msgtxt, "Found %d VMD MPI node%s containing a total of %d CPU%s and %d GPU%s:",
00132 numnodes, (numnodes > 1) ? "s" : "",
00133 totalcpus, (totalcpus > 1) ? "s" : "",
00134 totalgpus, (totalgpus != 1) ? "s" : "");
00135 msgInfo << msgtxt << sendmsg;
00136
00137 for (i=0; i<numnodes; i++) {
00138 sprintf(msgtxt,
00139 "%4d: %2d CPUs, %3.2fGB (%2ld%%) free mem, "
00140 "%d GPUs, "
00141
00142 "Name: %s",
00143 i, nodes[i].numcpus,
00144 nodes[i].corefree / 1024.0f, nodes[i].corepcnt,
00145 nodes[i].numgpus,
00146
00147 nodes[i].machname);
00148 msgInfo << msgtxt << sendmsg;;
00149 }
00150 }
00151
00152
00153
00154 MPI_Barrier(MPI_COMM_WORLD);
00155
00156 *noderank=mynode;
00157 *nodecount=numnodes;
00158
00159 return 0;
00160 }
00161
00162 #endif