예제 #1
0
void
MatlabExportModule :: doOutputHomogenizeDofIDs(TimeStep *tStep,    FILE *FID)
{

    std :: vector <FloatArray*> HomQuantities;
    double Vol = 0.0;

    // Initialize vector of arrays constaining homogenized quantities
    HomQuantities.resize(internalVarsToExport.giveSize());

    for (int j=0; j<internalVarsToExport.giveSize(); j++) {
        HomQuantities.at(j) = new FloatArray;
    }

    int nelem = this->elList.giveSize();
    for (int i = 1; i<=nelem; i++) {
        Element *e = this->emodel->giveDomain(1)->giveElement(elList.at(i));
        FEInterpolation *Interpolation = e->giveInterpolation();

        Vol = Vol + e->computeVolumeAreaOrLength();

        for ( GaussPoint *gp: *e->giveDefaultIntegrationRulePtr() ) {

            for (int j=0; j<internalVarsToExport.giveSize(); j++) {
                FloatArray elementValues;
                e->giveIPValue(elementValues, gp, (InternalStateType) internalVarsToExport(j), tStep);
                double detJ=fabs(Interpolation->giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e)));

                elementValues.times(gp->giveWeight()*detJ);
                if (HomQuantities.at(j)->giveSize() == 0) {
                    HomQuantities.at(j)->resize(elementValues.giveSize());
                    HomQuantities.at(j)->zero();
                };
                HomQuantities.at(j)->add(elementValues);
            }
        }
    }


    if (noscaling) Vol=1.0;

    for ( std :: size_t i = 0; i < HomQuantities.size(); i ++) {
        FloatArray *thisIS;
        thisIS = HomQuantities.at(i);
        thisIS->times(1.0/Vol);
        fprintf(FID, "\tspecials.%s = [", __InternalStateTypeToString ( InternalStateType (internalVarsToExport(i)) ) );

        for (int j = 0; j<thisIS->giveSize(); j++) {
            fprintf(FID, "%e", thisIS->at(j+1));
            if (j!=(thisIS->giveSize()-1) ) {
                fprintf(FID, ", ");
            }
        }
        fprintf(FID, "];\n");
        delete HomQuantities.at(i);
    }

}
예제 #2
0
void StructuralMaterialEvaluator :: doStepOutput(TimeStep *tStep)
{
    FloatArray outputValue;
    Domain *d = this->giveDomain(1);
    if ( tStep->isTheFirstStep() ) {
        this->outfile << "# Time";
        for ( int var : this->vars ) {
            this->outfile << ", " << __InternalStateTypeToString( ( InternalStateType ) var );
        }

        this->outfile << '\n';
    }

    outfile << tStep->giveIntrinsicTime();
    for ( int i = 1; i <= d->giveNumberOfMaterialModels(); i++ ) {
        Material *mat = d->giveMaterial(i);
        for ( int var : this->vars ) {
            mat->giveIPValue(outputValue, gps[i-1].get(), ( InternalStateType ) var, tStep);
            outfile << " " << outputValue;
        }
    }

    outfile << std :: endl;
}
void
NodalAveragingRecoveryModel :: initRegionMap(IntArray &regionMap, IntArray &regionValSize, InternalStateType type)
{
    int nregions = this->giveNumberOfVirtualRegions();
    int ielem, nelem = domain->giveNumberOfElements();
    int i, elemVR, regionsSkipped = 0;
    Element *element;
    NodalAveragingRecoveryModelInterface *interface;

    regionMap.resize(nregions);
    regionMap.zero();
    regionValSize.resize(nregions);
    regionValSize.zero();

    // loop over elements and check if implement interface
    for ( ielem = 1; ielem <= nelem; ielem++ ) {
        element = domain->giveElement(ielem);
        if ( !element-> isActivated(domain->giveEngngModel()->giveCurrentStep()) ) {  //skip inactivated elements
            continue;
        }
#ifdef __PARALLEL_MODE
        if ( element->giveParallelMode() != Element_local ) {
            continue;
        }

#endif

        if ( ( interface =  ( NodalAveragingRecoveryModelInterface * ) element->
                           giveInterface(NodalAveragingRecoveryModelInterfaceType) ) == NULL ) {
            // If an element doesn't implement the interface, it is ignored.

            //regionsSkipped = 1;
            //regionMap.at( element->giveRegionNumber() ) = 1;
            continue;
        } else {
            if ((elemVR = this->giveElementVirtualRegionNumber(ielem))) { // test if elemVR is nonzero
                if ( regionValSize.at(elemVR) ) {
                    if ( regionValSize.at(elemVR) !=
                        interface->NodalAveragingRecoveryMI_giveDofManRecordSize(type) ) {
                        // This indicates a size mis-match between different elements, no choice but to skip the region.
                        regionMap.at(elemVR) = 1;
                        regionsSkipped = 1;
                    }
                } else {
                    regionValSize.at(elemVR) = interface->
                                            NodalAveragingRecoveryMI_giveDofManRecordSize(type);
                }
            }
        }
    }

    if ( regionsSkipped ) {
        OOFEM_LOG_RELEVANT("NodalAveragingRecoveryModel :: initRegionMap: skipping regions for InternalStateType %s\n", __InternalStateTypeToString(type) );
        for ( i = 1; i <= nregions; i++ ) {
            if ( regionMap.at(i) ) {
                OOFEM_LOG_RELEVANT("%d ", i);
            }
        }

        OOFEM_LOG_RELEVANT("\n");
    }
}
예제 #4
0
int
SPRNodalRecoveryModel :: recoverValues(Set elementSet, InternalStateType type, TimeStep *tStep)
{
    int nnodes = domain->giveNumberOfDofManagers();
    IntArray regionNodalNumbers(nnodes);
    IntArray patchElems, dofManToDetermine, pap;
    FloatMatrix a;
    FloatArray dofManValues;
    IntArray dofManPatchCount;

    if ( ( this->valType == type ) && ( this->stateCounter == tStep->giveSolutionStateCounter() ) ) {
        return 1;
    }

#ifdef __PARALLEL_MODE
    this->initCommMaps();
#endif

    // clear nodal table
    this->clear();

    int regionValSize;
    int regionDofMans;


    regionValSize = 0;
    // loop over elements and determine local region node numbering and determine and check nodal values size
    if ( this->initRegionNodeNumbering(regionNodalNumbers, regionDofMans, elementSet) == 0 ) {
        return 0;
    }

    SPRPatchType regType = this->determinePatchType(elementSet);

    dofManPatchCount.resize(regionDofMans);
    dofManPatchCount.zero();


    //pap = patch assembly points
    this->determinePatchAssemblyPoints(pap, regType, elementSet);

    int npap = pap.giveSize();
    for ( int ipap = 1; ipap <= npap; ipap++ ) {
        int papNumber = pap.at(ipap);
        int oldSize = regionValSize;

        this->initPatch(patchElems, dofManToDetermine, pap, papNumber, elementSet);
        this->computePatch(a, patchElems, regionValSize, regType, type, tStep);
        if ( oldSize == 0 ) {
            dofManValues.resize(regionDofMans * regionValSize);
            dofManValues.zero();
        }
        this->determineValuesFromPatch(dofManValues, dofManPatchCount, regionNodalNumbers,
                                       dofManToDetermine, a, regType);
    }

#ifdef __PARALLEL_MODE
    this->exchangeDofManValues(dofManValues, dofManPatchCount, regionNodalNumbers, regionValSize);
#endif

    // average  recovered values of active region
    bool abortFlag = false;
    for ( int i = 1; i <= nnodes; i++ ) {
        if ( regionNodalNumbers.at(i) &&
            ( ( domain->giveDofManager(i)->giveParallelMode() == DofManager_local ) ||
             ( domain->giveDofManager(i)->giveParallelMode() == DofManager_shared ) ) ) {
            int eq = ( regionNodalNumbers.at(i) - 1 ) * regionValSize;
            if ( dofManPatchCount.at( regionNodalNumbers.at(i) ) ) {
                for ( int j = 1; j <= regionValSize; j++ ) {
                    dofManValues.at(eq + j) /= dofManPatchCount.at( regionNodalNumbers.at(i) );
                }
            } else {
                OOFEM_WARNING("values of %s in dofmanager %d undetermined", __InternalStateTypeToString(type), i);

                for ( int j = 1; j <= regionValSize; j++ ) {
                    dofManValues.at(eq + j) = 0.0;
                }
                //abortFlag = true;
            }
        }

        if ( abortFlag ) {
            abort();
        }

        // update recovered values
        this->updateRegionRecoveredValues(regionNodalNumbers, regionValSize, dofManValues);
    }

    this->valType = type;
    this->stateCounter = tStep->giveSolutionStateCounter();
    return 1;
}
int
NodalAveragingRecoveryModel :: recoverValues(Set elementSet, InternalStateType type, TimeStep *tStep)
{
    int nnodes = domain->giveNumberOfDofManagers();
    IntArray regionNodalNumbers(nnodes);
    IntArray regionDofMansConnectivity;
    FloatArray lhs, val;


    if ( ( this->valType == type ) && ( this->stateCounter == tStep->giveSolutionStateCounter() ) ) {
        return 1;
    }

#ifdef __PARALLEL_MODE
    bool parallel = this->domain->giveEngngModel()->isParallel();
    if ( parallel ) {
        this->initCommMaps();
    }
#endif

    // clear nodal table
    this->clear();

    int regionValSize = 0;
    int regionDofMans;

    // loop over elements and determine local region node numbering and determine and check nodal values size
    if ( this->initRegionNodeNumbering(regionNodalNumbers, regionDofMans, elementSet) == 0 ) {
        return 0;
    }

    regionDofMansConnectivity.resize(regionDofMans);
    regionDofMansConnectivity.zero();

    IntArray elements = elementSet.giveElementList();
    // assemble element contributions
    for ( int i = 1; i <= elements.giveSize(); i++ ) {
        int ielem = elements.at(i);
        NodalAveragingRecoveryModelInterface *interface;
        Element *element = domain->giveElement(ielem);

        if ( element->giveParallelMode() != Element_local ) {
            continue;
        }

        // If an element doesn't implement the interface, it is ignored.
        if ( ( interface = static_cast< NodalAveragingRecoveryModelInterface * >
                           ( element->giveInterface(NodalAveragingRecoveryModelInterfaceType) ) ) == NULL ) {
            //abort();
            continue;
        }

        int elemNodes = element->giveNumberOfDofManagers();
        // ask element contributions
        for ( int elementNode = 1; elementNode <= elemNodes; elementNode++ ) {
            int node = element->giveDofManager(elementNode)->giveNumber();
            interface->NodalAveragingRecoveryMI_computeNodalValue(val, elementNode, type, tStep);
            // if the element cannot evaluate this variable, it is ignored
            if ( val.giveSize() == 0 ) {
                continue;
            } else if ( regionValSize == 0 ) {
                regionValSize = val.giveSize();
                lhs.resize(regionDofMans * regionValSize);
                lhs.zero();
            } else if ( val.giveSize() != regionValSize ) {
                OOFEM_LOG_RELEVANT("NodalAveragingRecoveryModel :: size mismatch for InternalStateType %s, ignoring all elements that doesn't use the size %d\n", __InternalStateTypeToString(type), regionValSize);
                continue;
            }
            int eq = ( regionNodalNumbers.at(node) - 1 ) * regionValSize;
            for ( int j = 1; j <= regionValSize; j++ ) {
                lhs.at(eq + j) += val.at(j);
            }

            regionDofMansConnectivity.at( regionNodalNumbers.at(node) )++;
        }
    } // end assemble element contributions

#ifdef __PARALLEL_MODE
    if ( parallel ) {
        this->exchangeDofManValues(lhs, regionDofMansConnectivity, regionNodalNumbers, regionValSize);
    }
#endif

    // solve for recovered values of active region
    for ( int inode = 1; inode <= nnodes; inode++ ) {
        if ( regionNodalNumbers.at(inode) ) {
            int eq = ( regionNodalNumbers.at(inode) - 1 ) * regionValSize;
            for ( int i = 1; i <= regionValSize; i++ ) {
                if ( regionDofMansConnectivity.at( regionNodalNumbers.at(inode) ) > 0 ) {
                    lhs.at(eq + i) /= regionDofMansConnectivity.at( regionNodalNumbers.at(inode) );
                } else {
                    OOFEM_WARNING("values of dofmanager %d undetermined", inode);
                    lhs.at(eq + i) = 0.0;
                }
            }
        }
    }

    // update recovered values
    this->updateRegionRecoveredValues(regionNodalNumbers, regionValSize, lhs);

    this->valType = type;
    this->stateCounter = tStep->giveSolutionStateCounter();
    return 1;
}
예제 #6
0
//keyword "vars" in OOFEM input file
void
VTKXMLExportModule :: exportIntVarAs(InternalStateType valID, InternalStateValueType type,
                                     IntArray &mapG2L, IntArray &mapL2G, int regionDofMans, int ireg,
#ifdef __VTK_MODULE
                                     vtkSmartPointer<vtkUnstructuredGrid> &stream,
#else
                                     FILE *stream,
#endif
                                     TimeStep *tStep)
{
    Domain *d = emodel->giveDomain(1);
    int inode;
    int j, jsize, valSize;
    FloatArray iVal(3), t(9);
    const FloatArray *val;
    IntArray regionVarMap;

    this->giveSmoother();

#ifdef __VTK_MODULE
    vtkSmartPointer<vtkDoubleArray> intVarArray = vtkSmartPointer<vtkDoubleArray>::New();
    intVarArray->SetName(__InternalStateTypeToString(valID));
#endif

    if ( type == ISVT_SCALAR ) {
#ifdef __VTK_MODULE
        intVarArray->SetNumberOfComponents(1);
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" format=\"ascii\"> ", __InternalStateTypeToString(valID) );
#endif
    } else if ( type == ISVT_VECTOR ) {
#ifdef __VTK_MODULE
        intVarArray->SetNumberOfComponents(3);
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"3\" format=\"ascii\"> ",
                __InternalStateTypeToString(valID) );
#endif
    } else if ( ( type == ISVT_TENSOR_S3 ) || ( type == ISVT_TENSOR_S3E ) ) {
#ifdef __VTK_MODULE
        intVarArray->SetNumberOfComponents(9);
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"9\" format=\"ascii\"> ",
                __InternalStateTypeToString(valID) );
#endif
    } else if ( type == ISVT_TENSOR_G ) {
#ifdef __VTK_MODULE
        intVarArray->SetNumberOfComponents(9);
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"9\" format=\"ascii\"> ",
                __InternalStateTypeToString(valID) );
#endif
    } else {
        fprintf( stderr, "VTKXMLExportModule::exportIntVarAs: unsupported variable type %s\n", __InternalStateTypeToString(valID) );
    }
#ifdef __VTK_MODULE
    intVarArray->SetNumberOfTuples(regionDofMans);
#endif


    if ( !( ( valID == IST_DisplacementVector ) || ( valID == IST_MaterialInterfaceVal ) ) ) {
        this->smoother->recoverValues(valID, tStep);
        this->smoother->giveRegionRecordMap(regionVarMap, ireg, valID);
    }

    for ( inode = 1; inode <= regionDofMans; inode++ ) {
        if ( valID == IST_DisplacementVector ) { ///@todo Why does this code exists here? DisplacementVector isn't a internal variable. And if its really desired, then why the special treatment for just displacements?
            iVal.resize(3);
            val = & iVal;
            for ( j = 1; j <= 3; j++ ) {
                iVal.at(j) = d->giveNode( mapL2G.at(inode) )->giveUpdatedCoordinate(j, tStep, EID_MomentumBalance, 1.0) -
                             d->giveNode( mapL2G.at(inode) )->giveCoordinate(j);
            }
        } else if ( valID == IST_MaterialInterfaceVal ) {
            MaterialInterface *mi = emodel->giveMaterialInterface(1);
            if ( mi ) {
                iVal.resize(1);
                val = & iVal;
                iVal.at(1) = mi->giveNodalScalarRepresentation( mapL2G.at(inode) );
            }
        } else {
            int found = this->smoother->giveNodalVector(val, mapL2G.at(inode), ireg);
            if ( !found ) {
                iVal.resize( regionVarMap.giveSize() );
                iVal.zero();
                val = & iVal;
                //OOFEM_WARNING2("VTKXMLExportModule::exportIntVars: smoothing error: invalid data in node %d", inode);
            }
        }

        valSize = val->giveSize();
        if ( type == ISVT_SCALAR ) {
#ifdef __VTK_MODULE
            intVarArray->SetTuple1(inode-1, valSize ? val->at(1) : 0.0);
#else
            fprintf( stream, "%e ", valSize ? val->at(1) : 0.0 );
#endif
        } else if ( type == ISVT_VECTOR ) {
            jsize = min( 3, valSize );
            for ( j = 1; j <= jsize; j++ ) {
#ifdef __VTK_MODULE
                intVarArray->SetComponent(inode-1, j-1, val->at(j));
#else
                fprintf( stream, "%e ", val->at(j) );
#endif
            }

            for ( j = jsize + 1; j <= 3; j++ ) {
#ifdef __VTK_MODULE
                intVarArray->SetComponent(inode-1, j-1, 0.0);
#else
                fprintf(stream, "0.0 ");
#endif
            }
#ifndef __VTK_MODULE
            fprintf(stream, " ");
#endif

        } else if ( type == ISVT_TENSOR_S3 || type == ISVT_TENSOR_S3E ) {
            this->makeFullForm(t, *val, type, regionVarMap);

            for ( j = 1; j <= 9; j++ ) {
#ifdef __VTK_MODULE
                intVarArray->SetComponent(inode-1, j-1, t.at(j));
#else
                fprintf( stream, "%e ", t.at(j) );
#endif
            }
#ifndef __VTK_MODULE
            fprintf(stream, " ");
#endif
        } else if ( type == ISVT_TENSOR_G ) { // export general tensor values as scalars
            jsize = min( 9, val->giveSize() );
            for ( j = 1; j <= jsize; j++ ) {
#ifdef __VTK_MODULE
                intVarArray->SetComponent(inode-1, j-1, val->at(j));
#else
                fprintf( stream, "%e ", val->at(j) );
#endif

            }

            for ( j = jsize + 1; j <= 9; j++ ) {
#ifdef __VTK_MODULE
                intVarArray->SetComponent(inode-1, j-1, 0.0);
#else
                fprintf(stream, "0.0 ");
#endif
            }
#ifndef __VTK_MODULE
            fprintf(stream, " ");
#endif
        }
    } // end loop over dofmans

#ifdef __VTK_MODULE
    if (type == ISVT_SCALAR) {
        stream->GetPointData()->SetActiveScalars(__InternalStateTypeToString(valID));
        stream->GetPointData()->SetScalars(intVarArray);
    } else if (type == ISVT_VECTOR) {
        stream->GetPointData()->SetActiveVectors(__InternalStateTypeToString(valID));
        stream->GetPointData()->SetVectors(intVarArray);
    } else {
        stream->GetPointData()->SetActiveTensors(__InternalStateTypeToString(valID));
        stream->GetPointData()->SetTensors(intVarArray);
    }
#else
    fprintf(stream, "</DataArray>\n");
#endif
}
예제 #7
0
void
VTKXMLExportModule :: doOutput(TimeStep *tStep)
{
    if ( !testTimeStepOutput(tStep) ) {
        return;
    }
#ifdef __VTK_MODULE
    vtkSmartPointer<vtkUnstructuredGrid> stream = vtkSmartPointer<vtkUnstructuredGrid>::New();
    vtkSmartPointer<vtkPoints> nodes = vtkSmartPointer<vtkPoints>::New();
    vtkSmartPointer<vtkIdList> elemNodeArray = vtkSmartPointer<vtkIdList>::New();
#else
    FILE *stream = this->giveOutputStream(tStep);
    struct tm *current;
    time_t now;
    time(&now);
    current = localtime(&now);
    fprintf(stream, "<!-- TimeStep %e Computed %d-%02d-%02d at %02d:%02d:%02d -->\n", tStep->giveIntrinsicTime(), current->tm_year+1900, current->tm_mon+1, current->tm_mday, current->tm_hour,  current->tm_min,  current->tm_sec);
    fprintf(stream, "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");
    fprintf(stream, "<UnstructuredGrid>\n");
#endif


    Domain *d  = emodel->giveDomain(1);
    Element *elem;
    FloatArray *coords;
    int nelem = d->giveNumberOfElements();

    this->giveSmoother(); // make sure smoother is created

    // output nodes Region By Region
    int nregions = this->smoother->giveNumberOfVirtualRegions();
    int regionDofMans, totalcells;
    IntArray mapG2L, mapL2G;

    /* loop over regions */
    for ( int ireg = 1; ireg <= nregions; ireg++ ) {
        if ( ( ireg > 0 ) && ( this->regionsToSkip.contains(ireg) ) ) {
            continue;
        }

        // assemble local->global and global->local region map
        // and get number of single cells to process
        // the composite cells exported individually
        this->initRegionNodeNumbering(mapG2L, mapL2G, regionDofMans, totalcells, d, ireg);

        /* start default piece containing all single cell elements
         * the elements with composite geometry are assumed to be exported in individual pieces
         * after the default one
         */
#ifndef __PARALLEL_MODE
        if ( regionDofMans && totalcells ) {
#else
        if ( 1 ) {
#endif

#ifdef __VTK_MODULE
            for ( int inode = 1; inode <= regionDofMans; inode++ ) {
                coords = d->giveNode(mapL2G.at(inode))->giveCoordinates();
                int dims = coords->giveSize();
                nodes->InsertNextPoint(coords->at(1), dims >= 2 ? coords->at(2) : 0.0, dims >= 3 ? coords->at(3) : 0.0);
            }
            stream->SetPoints(nodes);
#else
            fprintf(stream, "<Piece NumberOfPoints=\"%d\" NumberOfCells=\"%d\">\n", regionDofMans, totalcells);
            // export nodes in region as vtk vertices
            fprintf(stream, "<Points>\n <DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\"> ");
            for ( int inode = 1; inode <= regionDofMans; inode++ ) {
                coords = d->giveNode( mapL2G.at(inode) )->giveCoordinates();
                for ( int i = 1; i <= coords->giveSize(); i++ ) {
                    fprintf( stream, "%e ", coords->at(i) );
                }

                for ( int i = coords->giveSize() + 1; i <= 3; i++ ) {
                    fprintf(stream, "%e ", 0.0);
                }
            }
            fprintf(stream, "</DataArray>\n</Points>\n");
#endif


            // output all cells of the piece
            int nelemNodes;
            IntArray cellNodes;
#ifdef __VTK_MODULE
            stream->Allocate(nelem);
#else
            fprintf(stream, "<Cells>\n");
            // output the connectivity data
            fprintf(stream, " <DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\"> ");
#endif
            for ( int ielem = 1; ielem <= nelem; ielem++ ) {
                elem = d->giveElement(ielem);
                if ( ( ireg > 0 ) && ( this->smoother->giveElementVirtualRegionNumber(ielem) != ireg ) ) {
                    continue;
                }

                if ( this->isElementComposite(elem) ) {
                    continue;                                  // composite cells exported individually
                }

                if ( !elem-> isActivated(tStep) ) {                  //skip inactivated elements
                    continue;
                }
                
#ifdef __PARALLEL_MODE
                if ( elem->giveParallelMode() != Element_local ) {
                    continue;
                }
#endif


                nelemNodes = elem->giveNumberOfNodes();
                this->giveElementCell(cellNodes, elem, 0);
#ifdef __VTK_MODULE
                elemNodeArray->Reset();
                elemNodeArray->SetNumberOfIds(nelemNodes);
#endif
                for ( int i = 1; i <= nelemNodes; i++ ) {
#ifdef __VTK_MODULE
                    elemNodeArray->SetId(i-1, mapG2L.at( cellNodes.at(i) ) - 1);
#else
                    fprintf(stream, "%d ", mapG2L.at( cellNodes.at(i) ) - 1);
#endif
                }
#ifdef __VTK_MODULE
                stream->InsertNextCell(this->giveCellType(elem), elemNodeArray);
#else
                fprintf(stream, " ");
#endif
            }

#ifndef __VTK_MODULE
            int vtkCellType;
            fprintf(stream, "</DataArray>\n");
            // output the offsets (index of individual element data in connectivity array)
            fprintf(stream, " <DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\"> ");
            int offset = 0;
            for ( int ielem = 1; ielem <= nelem; ielem++ ) {
                elem = d->giveElement(ielem);
                if ( ( ireg > 0 ) && ( this->smoother->giveElementVirtualRegionNumber(ielem) != ireg ) ) {
                    continue;
                }

#ifdef __PARALLEL_MODE
                if ( elem->giveParallelMode() != Element_local ) {
                    continue;
                }

#endif
                offset += elem->giveNumberOfNodes();
                fprintf(stream, "%d ", offset);
            }

            fprintf(stream, "</DataArray>\n");

            // output cell (element) types
            fprintf(stream, " <DataArray type=\"UInt8\" Name=\"types\" format=\"ascii\"> ");
            for ( int ielem = 1; ielem <= nelem; ielem++ ) {
                elem = d->giveElement(ielem);
                if ( ( ireg > 0 ) && ( this->smoother->giveElementVirtualRegionNumber(ielem) != ireg ) ) {
                    continue;
                }

                if ( this->isElementComposite(elem) ) {
                    continue;                                  // composite cells exported individually
                }

#ifdef __PARALLEL_MODE
                if ( elem->giveParallelMode() != Element_local ) {
                    continue;
                }

#endif
                vtkCellType = this->giveCellType(elem);
                fprintf(stream, "%d ", vtkCellType);
            }

            fprintf(stream, "</DataArray>\n");
            fprintf(stream, "</Cells>\n");
#endif

            // export primary and internal variables
#ifndef __VTK_MODULE
            this->exportPointDataHeader(stream, tStep);
#endif
            this->exportPrimaryVars(stream, mapG2L, mapL2G, regionDofMans, ireg, tStep);
            this->exportIntVars(stream, mapG2L, mapL2G, regionDofMans, ireg, tStep);
#ifndef __VTK_MODULE
            fprintf(stream, "</PointData>\n");
#endif

            //export cell data
            this->exportCellVars(stream, ireg, tStep);

#ifndef __VTK_MODULE
            // end of piece record
            fprintf(stream, "</Piece>\n");
#endif
        } // end of default piece for simple geometry elements

#if 1
        // loop over region elements with multi-cell geometry
        for ( int ielem = 1; ielem <= nelem; ielem++ ) {
            elem = d->giveElement(ielem);

            if ( this->regionsToSkip.contains( this->smoother->giveElementVirtualRegionNumber(ielem) ) ) {
                continue;
            }

 #ifdef __PARALLEL_MODE
            if ( elem->giveParallelMode() != Element_local ) {
                continue;
            }

 #endif

            if ( this->isElementComposite(elem) ) {
#ifndef __VTK_MODULE
                ///@todo Not sure how to deal with this.
                // multi cell (composite) elements should support vtkxmlexportmoduleinterface
                // and are exported as individual pieces (see VTKXMLExportModuleElementInterface)
                VTKXMLExportModuleElementInterface *interface =
                    ( VTKXMLExportModuleElementInterface * ) elem->giveInterface(VTKXMLExportModuleElementInterfaceType);
                if ( interface ) {
                    // passing this to access general piece related methods like exportPointDataHeader, etc.
                    interface->_export(stream, this, primaryVarsToExport, internalVarsToExport, tStep);
                }
#endif
            }
        } // end loop over multi-cell elements
#endif
    } // end loop over regions

    std::string fname = giveOutputFileName(tStep);

#ifdef __VTK_MODULE
#if 0
    // Doesn't as well as I would want it to, interface to VTK is to limited to control this.
    // * The PVTU-file is written by every process (seems to be impossible to avoid).
    // * Part files are renamed and time step and everything else is cut off => name collisions
    vtkSmartPointer<vtkXMLPUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
    writer->SetTimeStep(tStep->giveNumber()-1);
    writer->SetNumberOfPieces( this->emodel->giveNumberOfProcesses() );
    writer->SetStartPiece( this->emodel->giveRank() );
    writer->SetEndPiece( this->emodel->giveRank() );
#else
    vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
#endif
    writer->SetFileName(fname.c_str());
    writer->SetInput(stream);
    // Optional - set the mode. The default is binary.
    //writer->SetDataModeToBinary();
    //writer->SetDataModeToAscii();
    writer->Write();
#else
    // finish unstructured grid data and vtk file
    fprintf(stream, "</UnstructuredGrid>\n</VTKFile>");
    fclose(stream);
#endif

    // Writing the *.pvd-file. Only time step information for now. It's named "timestep" but is actually the total time.
    // First we check to see that there are more than 1 time steps, otherwise it is redundant;
#ifdef __PARALLEL_MODE
    if ( emodel->isParallel() && emodel->giveRank() == 0 ) {
        ///@todo Should use probably use PVTU-files instead. It is starting to get messy.
        // For this to work, all processes must have an identical output file name.
        for (int i = 0; i < this->emodel->giveNumberOfProcesses(); ++i) {
            std::ostringstream pvdEntry;
            char fext[100];
            if (this->emodel->giveNumberOfProcesses() > 1) {
                sprintf( fext, "_%03d.m%d.%d", i, this->number, tStep->giveNumber() );
            } else {
                sprintf( fext, "m%d.%d", this->number, tStep->giveNumber() );
            }
            pvdEntry << "<DataSet timestep=\"" << tStep->giveIntrinsicTime() << "\" group=\"\" part=\"" << i << "\" file=\""
                    << this->emodel->giveOutputBaseFileName() << fext << ".vtu\"/>";
            this->pvdBuffer.push_back(pvdEntry.str());
        }
        this->writeVTKCollection();
    } else
#endif
    if ( !emodel->isParallel() && tStep->giveNumber() >= 1 ) { // For non-parallel enabled OOFEM, then we only check for multiple steps.
        std::ostringstream pvdEntry;
        pvdEntry << "<DataSet timestep=\"" << tStep->giveIntrinsicTime() << "\" group=\"\" part=\"\" file=\"" << fname << "\"/>";
        this->pvdBuffer.push_back(pvdEntry.str());
        this->writeVTKCollection();
    }
}

#ifndef __VTK_MODULE
void
VTKXMLExportModule :: exportPointDataHeader(FILE *stream, TimeStep *tStep)
{
    int n;
    std :: string scalars, vectors, tensors;

    n = primaryVarsToExport.giveSize();

    UnknownType type;

    for ( int i = 1; i <= n; i++ ) {
        type = ( UnknownType ) primaryVarsToExport.at(i);
        if ( ( type == DisplacementVector ) || ( type == EigenVector ) || ( type == VelocityVector ) ) {
            vectors += __UnknownTypeToString(type);
            vectors.append(" ");
        } else if ( ( type == FluxVector ) || ( type == PressureVector ) || ( type == Temperature ) ) {
            scalars += __UnknownTypeToString(type);
            scalars.append(" ");
        } else {
            OOFEM_ERROR2( "VTKXMLExportModule::exportPrimVarAs: unsupported UnknownType %s", __UnknownTypeToString(type) );
        }
    }

    InternalStateType isttype;
    InternalStateValueType vtype;


    n = internalVarsToExport.giveSize();

    // prepare header
    for ( int i = 1; i <= n; i++ ) {
        isttype = ( InternalStateType ) internalVarsToExport.at(i);
        vtype = giveInternalStateValueType(isttype);

        if ( vtype == ISVT_SCALAR ) {
            scalars += __InternalStateTypeToString(isttype);
            scalars.append(" ");
        } else if ( vtype == ISVT_VECTOR ) {
            vectors += __InternalStateTypeToString(isttype);
            vectors.append(" ");
        } else if ( ( vtype == ISVT_TENSOR_S3 ) || ( vtype == ISVT_TENSOR_S3E ) ) {
            tensors += __InternalStateTypeToString(isttype);
            tensors.append(" ");
        } else if ( vtype == ISVT_TENSOR_G ) {
            vectors += __InternalStateTypeToString(isttype);
            vectors.append(" ");
        } else {
            fprintf( stderr, "VTKXMLExportModule::exportIntVars: unsupported variable type %s\n", __InternalStateTypeToString(isttype) );
        }
    }

    // print header
    fprintf( stream, "<PointData Scalars=\"%s\" Vectors=\"%s\" Tensors=\"%s\" >\n",
            scalars.c_str(), vectors.c_str(), tensors.c_str() );
}
예제 #8
0
//keyword "cellvars" in OOFEM input file
void
VTKXMLExportModule :: exportCellVarAs(InternalStateType type, int region,
#ifdef __VTK_MODULE
    vtkSmartPointer<vtkUnstructuredGrid> &stream,
#else
    FILE *stream,
#endif
    TimeStep *tStep)
{
    Domain *d = emodel->giveDomain(1);
    int ielem, nelem = d->giveNumberOfElements();
    int pos;
    Element *elem;
    FloatMatrix mtrx(3, 3);
    IntegrationRule *iRule;
    GaussPoint *gp;
    FloatArray answer, temp;
    double gptot;
    int ncomponents = 1;

#ifdef __VTK_MODULE
    vtkSmartPointer<vtkDoubleArray> cellVarsArray = vtkSmartPointer<vtkDoubleArray>::New();
    cellVarsArray->SetName(__InternalStateTypeToString(type));
#endif

    switch ( type ) {
    case IST_MaterialNumber:
    case IST_ElementNumber:
    case IST_Pressure:
        // if it wasn't for IST_Pressure,
#ifdef __VTK_MODULE
        cellVarsArray->SetNumberOfComponents(1);
        cellVarsArray->SetNumberOfTuples(nelem);
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" format=\"ascii\">\n", __InternalStateTypeToString(type) );
#endif
        for ( ielem = 1; ielem <= nelem; ielem++ ) {
            elem = d->giveElement(ielem);

            if ( (( region > 0 ) && ( this->smoother->giveElementVirtualRegionNumber(ielem) != region ))
                    || this->isElementComposite(elem) || !elem-> isActivated(tStep) ) { // composite cells exported individually
                continue;
            }

#ifdef __PARALLEL_MODE
            if ( elem->giveParallelMode() != Element_local ) {
                continue;
            }

#endif
            if ( type == IST_MaterialNumber ) {
#ifdef __VTK_MODULE
                cellVarsArray->SetTuple1(ielem-1, elem->giveMaterial()->giveNumber() ); // Should be integer..
#else
                fprintf( stream, "%d ", elem->giveMaterial()->giveNumber() );
#endif
            } else if ( type == IST_ElementNumber ) {
#ifdef __VTK_MODULE
                cellVarsArray->SetTuple1(ielem-1,  elem->giveNumber() ); // Should be integer..
#else
                fprintf( stream, "%d ", elem->giveNumber() );
#endif
            } else if (type == IST_Pressure) { ///@todo Why this special treatment for pressure? / Mikael
                if (elem->giveNumberOfInternalDofManagers() == 1) {
                    IntArray pmask(1); pmask.at(1) = P_f;
                    elem->giveInternalDofManager(1)->giveUnknownVector (answer, pmask,EID_ConservationEquation, VM_Total, tStep);
#ifdef __VTK_MODULE
                    cellVarsArray->SetTuple1(ielem-1,  answer.at(1) ); // Should be integer..
#else
                    fprintf( stream, "%f ", answer.at(1) );
#endif
                }
            }
        }
#ifdef __VTK_MODULE
        stream->GetCellData()->SetActiveScalars(__InternalStateTypeToString(type));
        stream->GetCellData()->SetScalars(cellVarsArray);
#endif
        break;

    case IST_MaterialOrientation_x:
    case IST_MaterialOrientation_y:
    case IST_MaterialOrientation_z:
#ifdef __VTK_MODULE
        cellVarsArray->SetNumberOfComponents(3);
        cellVarsArray->SetNumberOfTuples(nelem);
        ncomponents = 3;
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"3\" format=\"ascii\">\n", __InternalStateTypeToString(type) );
#endif
        if ( type == IST_MaterialOrientation_x ) {
            pos = 1;
        }

        if ( type == IST_MaterialOrientation_y ) {
            pos = 2;
        }

        if ( type == IST_MaterialOrientation_z ) {
            pos = 3;
        }

        for ( ielem = 1; ielem <= nelem; ielem++ ) {
            ///@todo Should no elements be skipped here? / Mikael
            if ( !d->giveElement(ielem)->giveLocalCoordinateSystem(mtrx) ) {
                mtrx.resize(3, 3);
                mtrx.zero();
            }
#ifdef __VTK_MODULE
            cellVarsArray->SetTuple3(ielem-1,  mtrx.at(1, pos), mtrx.at(2,pos), mtrx.at(3,pos) );
#else
            fprintf( stream, "%f %f %f  ", mtrx.at(1, pos), mtrx.at(2, pos), mtrx.at(3, pos) );
#endif
        }

#ifdef __VTK_MODULE
        stream->GetCellData()->SetActiveVectors(__InternalStateTypeToString(type));
        stream->GetCellData()->SetVectors(cellVarsArray);
#endif
        break;

    default:
        bool reshape = false;
        InternalStateValueType vt = giveInternalStateValueType(type);
        if ( vt == ISVT_SCALAR ) {
            ncomponents = 1;
        } else if ( vt == ISVT_VECTOR ) {
            ncomponents = 3;
        } else {
            ncomponents = 9;
            reshape = true;
        }
#ifdef __VTK_MODULE
        cellVarsArray->SetNumberOfComponents(ncomponents);
        cellVarsArray->SetNumberOfTuples(nelem);
#else
        fprintf( stream, "<DataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"%d\" format=\"ascii\">\n", __InternalStateTypeToString(type), ncomponents );
#endif

        IntArray redIndx;
        for ( ielem = 1; ielem <= nelem; ielem++ ) {
            elem = d->giveElement(ielem);
            if ( (( region > 0 ) && ( this->smoother->giveElementVirtualRegionNumber(ielem) != region ))
                    || this->isElementComposite(elem) || !elem-> isActivated(tStep) ) { // composite cells exported individually
                continue;
            }
#ifdef __PARALLEL_MODE
            if ( elem->giveParallelMode() != Element_local ) {
                continue;
            }
#endif
            gptot = 0;
            answer.resize(0);
            iRule = elem->giveDefaultIntegrationRulePtr();
            if (iRule) {
                MaterialMode mmode = _Unknown;
                for (int i = 0; i < iRule->getNumberOfIntegrationPoints(); ++i) {
                    gp = iRule->getIntegrationPoint(i);
                    mmode = gp->giveMaterialMode();
                    elem->giveIPValue(temp, gp, type, tStep);
                    gptot += gp->giveWeight();
                    answer.add(gp->giveWeight(), temp);
                }
                answer.times(1./gptot);
                elem->giveMaterial()->giveIntVarCompFullIndx(redIndx, type, mmode);
            }
            // Reshape the Voigt vectors to include all components (duplicated if necessary, VTK insists on 9 components for tensors.)
            if ( reshape && answer.giveSize() != 9) { // If it has 9 components, then it is assumed to be proper already.
                FloatArray tmp = answer;
                this->makeFullForm(answer, tmp, vt, redIndx);
            } else if ( vt == ISVT_VECTOR && answer.giveSize() < 3) {
                answer.setValues(3,
                                 answer.giveSize() > 1 ? answer.at(1) : 0.0,
                                 answer.giveSize() > 2 ? answer.at(2) : 0.0,
                                 0.0);
            } else if ( ncomponents != answer.giveSize() ) { // Trying to gracefully handle bad cases, just output zeros.
                answer.resize(ncomponents);
                answer.zero();
            }
            for (int i = 1; i <= ncomponents; ++i) {
#ifdef __VTK_MODULE
                cellVarsArray->SetComponent(ielem-1, i-1, answer.at(i));
#else
                fprintf( stream, "%e ", answer.at(i) );
#endif
            }
#ifndef __VTK_MODULE
            fprintf( stream, "\n" );
#endif
        }
#ifdef __VTK_MODULE
        if (ncomponents == 1) {
            stream->GetCellData()->SetActiveScalars(__InternalStateTypeToString(type));
            stream->GetCellData()->SetScalars(cellVarsArray);
        } else if (ncomponents == 3) {
            stream->GetCellData()->SetActiveVectors(__InternalStateTypeToString(type));
            stream->GetCellData()->SetVectors(cellVarsArray);
        } else {
            stream->GetCellData()->SetActiveTensors(__InternalStateTypeToString(type));
            stream->GetCellData()->SetTensors(cellVarsArray);
        }
#endif
    }
#ifndef __VTK_MODULE
    fprintf(stream, "</DataArray>\n");
#endif
}