bool Foam::movingBodyTopoFvMesh::update()
{
    // Store points to recreate mesh motion
    pointField oldPointsNew = allPoints();
    pointField newPoints = allPoints();

    autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();

    bool localMeshChanged = topoChangeMap->morphing();
    bool globalMeshChanged = localMeshChanged;
    reduce(globalMeshChanged, orOp<bool>());

    if (globalMeshChanged)
    {
        Pout<< "Topology change. Calculating motion point mask" << endl;
        motionMask_ = calcMotionMask();
    }

    if (localMeshChanged)
    {
//         // Map old points onto the new mesh
//         pointField mappedOldPointsNew(allPoints().size());
//         mappedOldPointsNew.map(oldPointsNew, topoChangeMap->pointMap());

//         movePoints(mappedOldPointsNew);
//         resetMotion();
//         setV0();

        // Get new points from preMotion
        newPoints = topoChangeMap().preMotionPoints();
    }
//     else
//     {
//         // No change, use old points
//         movePoints(oldPointsNew);
//         resetMotion();
//         setV0();
//     }

    // Calculate new points using a velocity transformation
    newPoints += motionMask_*
        transform(SBMFPtr_().velocity(), newPoints)*time().deltaT().value();

    Info << "Executing mesh motion" << endl;
    movePoints(newPoints);

    return localMeshChanged;
}
Пример #2
0
bool Foam::rawTopoChangerFvMesh::update()
{
    // Do mesh changes (use inflation - put new points in topoChangeMap)
    Info<< "rawTopoChangerFvMesh : Checking for topology changes..."
        << endl;

    // Mesh not moved/changed yet
    moving(false);
    topoChanging(false);

    // Do any topology changes. Sets topoChanging (through polyTopoChange)
    autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);

    bool hasChanged = topoChangeMap.valid();

    if (hasChanged)
    {
        Info<< "rawTopoChangerFvMesh : Done topology changes..."
            << endl;

        // Temporary: fix fields on patch faces created out of nothing
        // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

        // Two situations:
        // - internal faces inflated out of nothing
        // - patch faces created out of previously internal faces

        // Is face mapped in any way?
        PackedBoolList mappedFace(nFaces());

        const label nOldInternal = topoChangeMap().oldPatchStarts()[0];

        const labelList& faceMap = topoChangeMap().faceMap();
        for (label facei = 0; facei < nInternalFaces(); facei++)
        {
            if (faceMap[facei] >= 0)
            {
                mappedFace[facei] = 1;
            }
        }
        for (label facei = nInternalFaces(); facei < nFaces(); facei++)
        {
            if (faceMap[facei] >= 0 && faceMap[facei] >= nOldInternal)
            {
                mappedFace[facei] = 1;
            }
        }

        const List<objectMap>& fromFaces = topoChangeMap().facesFromFacesMap();

        forAll(fromFaces, i)
        {
            mappedFace[fromFaces[i].index()] = 1;
        }

        const List<objectMap>& fromEdges = topoChangeMap().facesFromEdgesMap();

        forAll(fromEdges, i)
        {
            mappedFace[fromEdges[i].index()] = 1;
        }

        const List<objectMap>& fromPts = topoChangeMap().facesFromPointsMap();

        forAll(fromPts, i)
        {
            mappedFace[fromPts[i].index()] = 1;
        }

        // Set unmapped faces to zero
        Info<< "rawTopoChangerFvMesh : zeroing unmapped boundary values."
            << endl;
        zeroUnmappedValues<scalar, fvPatchField, volMesh>(mappedFace);
        zeroUnmappedValues<vector, fvPatchField, volMesh>(mappedFace);
        zeroUnmappedValues<sphericalTensor, fvPatchField, volMesh>(mappedFace);
        zeroUnmappedValues<symmTensor, fvPatchField, volMesh>(mappedFace);
        zeroUnmappedValues<tensor, fvPatchField, volMesh>(mappedFace);

        // Special handling for phi: set unmapped faces to recreated phi
        Info<< "rawTopoChangerFvMesh :"
            << " recreating phi for unmapped boundary values." << endl;
        const volVectorField& U = lookupObject<volVectorField>("U");
        surfaceScalarField& phi = const_cast<surfaceScalarField&>
        (
            lookupObject<surfaceScalarField>("phi")
        );
        setUnmappedValues
        (
            phi,
            mappedFace,
            (linearInterpolate(U) & Sf())()
        );


        if (topoChangeMap().hasMotionPoints())
        {
            pointField newPoints = topoChangeMap().preMotionPoints();

            // Give the meshModifiers opportunity to modify points
            Info<< "rawTopoChangerFvMesh :"
                << " calling modifyMotionPoints." << endl;
            topoChanger_.modifyMotionPoints(newPoints);

            // Actually move points
            Info<< "rawTopoChangerFvMesh :"
                << " calling movePoints." << endl;

            movePoints(newPoints);
        }
    }
    else
    {
        //Pout<< "rawTopoChangerFvMesh :"
        //    << " no topology changes..." << endl;
    }

    return hasChanged;
Пример #3
0
bool Foam::movingConeTopoFvMesh::update()
{
    // Do mesh changes (use inflation - put new points in topoChangeMap)
    autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);

    // Calculate the new point positions depending on whether the
    // topological change has happened or not
    pointField newPoints;

    vector curMotionVel_ =
        motionVelAmplitude_*
        Foam::sin(time().value()*M_PI/motionVelPeriod_);

    Pout<< "time:" << time().value() << " curMotionVel_:" << curMotionVel_
        << " curLeft:" << curLeft_ << " curRight:" << curRight_
        << endl;

    if (topoChangeMap.valid())
    {
        Info<< "Topology change. Calculating motion points" << endl;

        if (topoChangeMap().hasMotionPoints())
        {
            Info<< "Topology change. Has premotion points" << endl;
            //Info<< "preMotionPoints:" << topoChangeMap().preMotionPoints()
            //    << endl;

            //mkDir(time().timePath());
            //{
            //    OFstream str(time().timePath()/"meshPoints.obj");
            //    Pout<< "Writing mesh with meshPoints to " << str.name()
            //        << endl;
            //
            //    const pointField& currentPoints = points();
            //    label vertI = 0;
            //    forAll(currentPoints, pointI)
            //    {
            //        meshTools::writeOBJ(str, currentPoints[pointI]);
            //        vertI++;
            //    }
            //    forAll(edges(), edgeI)
            //    {
            //        const edge& e = edges()[edgeI];
            //        str << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
            //    }
            //}
            //{
            //    OFstream str(time().timePath()/"preMotionPoints.obj");
            //    Pout<< "Writing mesh with preMotionPoints to " << str.name()
            //        << endl;
            //
            //    const pointField& newPoints =
            //        topoChangeMap().preMotionPoints();
            //    label vertI = 0;
            //    forAll(newPoints, pointI)
            //    {
            //        meshTools::writeOBJ(str, newPoints[pointI]);
            //        vertI++;
            //    }
            //    forAll(edges(), edgeI)
            //    {
            //        const edge& e = edges()[edgeI];
            //        str << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
            //    }
            //}


            motionMask_ =
                vertexMarkup
                (
                    topoChangeMap().preMotionPoints(),
                    curLeft_,
                    curRight_
                );

            // Move points inside the motionMask
            newPoints =
                topoChangeMap().preMotionPoints()
              + (
                    pos(0.5 - mag(motionMask_)) // cells above the body
                )*curMotionVel_*time().deltaT().value();
        }
        else
        {
            Info<< "Topology change. Already set mesh points" << endl;

            motionMask_ =
                vertexMarkup
                (
                    points(),
                    curLeft_,
                    curRight_
                );

            // Move points inside the motionMask
            newPoints =
                points()
              + (
                    pos(0.5 - mag(motionMask_)) // cells above the body
                )*curMotionVel_*time().deltaT().value();
        }
    }
    else
    {
        Info<< "No topology change" << endl;
        // Set the mesh motion
        newPoints =
            points()
          + (
                pos(0.5 - mag(motionMask_)) // cells above the body
           )*curMotionVel_*time().deltaT().value();
    }

    // The mesh now contains the cells with zero volume
    Info << "Executing mesh motion" << endl;
    movePoints(newPoints);
    //  The mesh now has got non-zero volume cells

    curLeft_ = average
    (
        faceZones()
        [
            faceZones().findZoneID("leftExtrusionFaces")
        ]().localPoints()
    ).x() - SMALL;

    curRight_ = average
    (
        faceZones()
        [
            faceZones().findZoneID("rightExtrusionFaces")
        ]().localPoints()
    ).x() + SMALL;


    return true;
}