void Spherical_MaxVel_Output( UnderworldContext* context ) {
   Spherical_CubedSphereNusselt* 	self		= (Spherical_CubedSphereNusselt*)LiveComponentRegister_Get( context->CF->LCRegister, (Name)Spherical_CubedSphereNusselt_Type );
   FeVariable*				velocityField   = self->velocityField;
   FeMesh* 				mesh 		= velocityField->feMesh;
   Grid*				grid		= NULL;
   unsigned*				sizes		= NULL;
   int 					vertId, dId;
   unsigned 				dT_i, dT_j, nDomainSize, ijk[3];
   double 				vel[3], velMax[2], gVelMax[2], velMag;

   //set 'em big and negative
   velMax[0] = velMax[1] = -1*HUGE_VAL;

   // get vert grid
   RegularMeshUtils_ErrorCheckAndGetDetails( (Mesh*)mesh, MT_VERTEX, &nDomainSize, &grid );
   sizes = grid->sizes;

   // go around 
   for( dT_i = 0; dT_i < sizes[1]; dT_i++ ) {
      for( dT_j = 0; dT_j < sizes[2]; dT_j++ ) {
         // find inner vertex
         ijk[0] = 0;
         // angular discretisation
         ijk[1] = dT_i;
         ijk[2] = dT_j;
         vertId = Grid_Project( grid, ijk );

         // if the node is local
         if( Mesh_GlobalToDomain( mesh, MT_VERTEX, vertId, &dId ) == True ) {
            FeVariable_GetValueAtNode( velocityField, dId, vel );
            velMag = sqrt( vel[0]*vel[0] + vel[1]*vel[1] + vel[2]*vel[2] );
            if( velMag > velMax[0] ) velMax[0] = velMag;
         }

         // find outer vertex
         ijk[0] = grid->sizes[0]-1;
         vertId = Grid_Project( grid, ijk );

         // if the node is local
         if( Mesh_GlobalToDomain( mesh, MT_VERTEX, vertId, &dId ) == True ) {
            FeVariable_GetValueAtNode( velocityField, dId, vel );
            velMag = sqrt( vel[0]*vel[0] + vel[1]*vel[1] + vel[2]*vel[2] );
            if( velMag > velMax[1] ) velMax[1] = velMag;
         }
      }
   }

   (void)MPI_Allreduce( velMax, gVelMax, 2, MPI_DOUBLE, MPI_MAX, context->communicator );

   StgFEM_FrequentOutput_PrintValue( context, gVelMax[1] ); // print outer velocity max
   StgFEM_FrequentOutput_PrintValue( context, gVelMax[0] ); // print inner velocity max
}
Ejemplo n.º 2
0
void lucMeshCrossSection_Sample( void* drawingObject, Bool reverse)
{
   lucMeshCrossSection* self          = (lucMeshCrossSection*)drawingObject;
   FeVariable*          fieldVariable = (FeVariable*) self->fieldVariable;
   Mesh*                mesh          = (Mesh*) fieldVariable->feMesh;
   Grid*                vertGrid;
   Node_LocalIndex      crossSection_I;
   IJK                  node_ijk;
   Node_GlobalIndex     node_gI;
   Node_DomainIndex     node_dI;
   int                  i,j, d, sizes[3] = {1,1,1};
   Coord                globalMin, globalMax, min, max;

      int localcount = 0;

   vertGrid = *(Grid**)ExtensionManager_Get( mesh->info, mesh, self->vertexGridHandle );
   for (d=0; d<fieldVariable->dim; d++) sizes[d] = vertGrid->sizes[d];
   self->dim[0] = sizes[ self->axis ];
   self->dim[1] = sizes[ self->axis1 ];
   self->dim[2] = sizes[ self->axis2 ];

   crossSection_I = lucCrossSection_GetValue(self, 0, self->dim[0]-1);

   FieldVariable_GetMinAndMaxLocalCoords( fieldVariable, min, max );
   FieldVariable_GetMinAndMaxGlobalCoords( fieldVariable, globalMin, globalMax );

   Journal_Printf( lucDebug, "%s called on field %s, with axis of cross section as %d, crossSection_I as %d (dims %d,%d,%d) field dim %d\n",
                    __func__, fieldVariable->name, self->axis, crossSection_I, self->dim[0], self->dim[1], self->dim[2], self->fieldDim);

   /* Get mesh cross section self->vertices and values */
   self->resolutionA = self->dim[1];
   self->resolutionB = self->dim[2];
   lucCrossSection_AllocateSampleData(self, self->fieldDim);
   int lSize = Mesh_GetLocalSize( mesh, MT_VERTEX );
   double time = MPI_Wtime();
   Journal_Printf(lucInfo, "Sampling mesh (%s) %d x %d...  0%", self->name, self->dim[1], self->dim[2]);
   node_ijk[ self->axis ] = crossSection_I;
   for ( i = 0 ; i < self->dim[1]; i++ )
   {
      int percent = 100 * (i + 1) / self->dim[1];
      Journal_Printf(lucInfo, "\b\b\b\b%3d%%", percent);
      fflush(stdout);

      /* Reverse order if requested */
      int i0 = i;
      if (reverse) i0 = self->dim[1] - i - 1;

      node_ijk[ self->axis1 ] = i0;

      for ( j = 0 ; j < self->dim[2]; j++ )
      {
         self->vertices[i][j][0] = HUGE_VAL;
         self->vertices[i][j][2] = 0;
         node_ijk[ self->axis2 ] = j;
         node_gI = Grid_Project( vertGrid, node_ijk );
         /* Get coord and value if node is local... */
         if (Mesh_GlobalToDomain( mesh, MT_VERTEX, node_gI, &node_dI ) && node_dI < lSize)
         {  
            /* Found on this processor */
            double value[self->fieldDim];
            FeVariable_GetValueAtNode( fieldVariable, node_dI, value );
            double* pos = Mesh_GetVertex( mesh, node_dI );
            /*fprintf(stderr, "[%d] (%d,%d) Node %d %f,%f,%f value %f\n", self->context->rank, i, j, node_gI, pos[0], pos[1], pos[2], value);*/
         
            for (d=0; d<fieldVariable->dim; d++)
               self->vertices[i][j][d] = pos[d];

            for (d=0; d<self->fieldDim; d++)
               self->values[i][j][d] = (float)value[d];

            localcount++;
         }
      }
   }
   Journal_Printf(lucInfo, " %f sec. ", MPI_Wtime() - time);

   /* Collate */
   time = MPI_Wtime();
   for ( i=0 ; i < self->dim[1]; i++ )
   {
      for ( j=0 ; j < self->dim[2]; j++ )
      {
         /* Receive values at root */
         if (self->context->rank == 0)
         {
            /* Already have value? */
            if (self->vertices[i][j][0] != HUGE_VAL) {localcount--; continue; }

            /* Recv (pos and value together = (3 + fevar dims)*float) */
            float data[3 + self->fieldDim];
            (void)MPI_Recv(data, 3+self->fieldDim, MPI_FLOAT, MPI_ANY_SOURCE, i*self->dim[2]+j, self->context->communicator, MPI_STATUS_IGNORE);
            /* Copy */
            memcpy(self->vertices[i][j], data, 3 * sizeof(float));
            memcpy(self->values[i][j], &data[3], self->fieldDim * sizeof(float));
         }
         else
         {
            /* Found on this proc? */
            if (self->vertices[i][j][0] == HUGE_VAL) continue;

            /* Copy */
            float data[3 + self->fieldDim];
            memcpy(data, self->vertices[i][j], 3 * sizeof(float));
            memcpy(&data[3], self->values[i][j], self->fieldDim * sizeof(float));

            /* Send values to root (pos & value = 4 * float) */
            MPI_Ssend(data, 3+self->fieldDim, MPI_FLOAT, 0, i*self->dim[2]+j, self->context->communicator);
            localcount--;
         }
      }
   }
   MPI_Barrier(self->context->communicator);    /* Barrier required, prevent subsequent MPI calls from interfering with transfer */
   Journal_Printf(lucInfo, " Gather in %f sec.\n", MPI_Wtime() - time);
   Journal_Firewall(localcount == 0, lucError,
                     "Error - in %s: count of values sampled compared to sent/received by mpi on proc %d does not match (balance = %d)\n",
                     __func__, self->context->rank, localcount);
}