writeFields.C 12.7 KB
Newer Older
1 2
#include "writeFields.H"
#include "volFields.H"
3
#include "surfaceFields.H"
4 5 6
#include "polyMeshTools.H"
#include "zeroGradientFvPatchFields.H"
#include "syncTools.H"
7
#include "tetPointRef.H"
8
#include "regionSplit.H"
9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71

using namespace Foam;

void maxFaceToCell
(
    const scalarField& faceData,
    volScalarField& cellData
)
{
    const cellList& cells = cellData.mesh().cells();

    scalarField& cellFld = cellData.ref();

    cellFld = -GREAT;
    forAll(cells, cellI)
    {
        const cell& cFaces = cells[cellI];
        forAll(cFaces, i)
        {
            cellFld[cellI] = max(cellFld[cellI], faceData[cFaces[i]]);
        }
    }

    forAll(cellData.boundaryField(), patchI)
    {
        fvPatchScalarField& fvp = cellData.boundaryFieldRef()[patchI];

        fvp = fvp.patch().patchSlice(faceData);
    }
    cellData.correctBoundaryConditions();
}


void minFaceToCell
(
    const scalarField& faceData,
    volScalarField& cellData
)
{
    const cellList& cells = cellData.mesh().cells();

    scalarField& cellFld = cellData.ref();

    cellFld = GREAT;
    forAll(cells, cellI)
    {
        const cell& cFaces = cells[cellI];
        forAll(cFaces, i)
        {
            cellFld[cellI] = min(cellFld[cellI], faceData[cFaces[i]]);
        }
    }

    forAll(cellData.boundaryField(), patchI)
    {
        fvPatchScalarField& fvp = cellData.boundaryFieldRef()[patchI];

        fvp = fvp.patch().patchSlice(faceData);
    }
    cellData.correctBoundaryConditions();
}


72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117
void minFaceToCell
(
    const surfaceScalarField& faceData,
    volScalarField& cellData,
    const bool correctBoundaryConditions
)
{
    scalarField& cellFld = cellData.ref();

    cellFld = GREAT;

    const labelUList& own = cellData.mesh().owner();
    const labelUList& nei = cellData.mesh().neighbour();

    // Internal faces
    forAll(own, facei)
    {
        cellFld[own[facei]] = min(cellFld[own[facei]], faceData[facei]);
        cellFld[nei[facei]] = min(cellFld[nei[facei]], faceData[facei]);
    }

    // Patch faces
    forAll(faceData.boundaryField(), patchi)
    {
        const fvsPatchScalarField& fvp = faceData.boundaryField()[patchi];
        const labelUList& fc = fvp.patch().faceCells();

        forAll(fc, i)
        {
            cellFld[fc[i]] = min(cellFld[fc[i]], fvp[i]);
        }
    }

    volScalarField::Boundary& bfld = cellData.boundaryFieldRef();

    forAll(bfld, patchi)
    {
        bfld[patchi] = faceData.boundaryField()[patchi];
    }
    if (correctBoundaryConditions)
    {
        cellData.correctBoundaryConditions();
    }
}


118 119 120
void Foam::writeFields
(
    const fvMesh& mesh,
121
    const wordHashSet& selectedFields
122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148
)
{
    if (selectedFields.empty())
    {
        return;
    }

    Info<< "Writing fields with mesh quality parameters" << endl;

    if (selectedFields.found("nonOrthoAngle"))
    {
        //- Face based orthogonality
        const scalarField faceOrthogonality
        (
            polyMeshTools::faceOrthogonality
            (
                mesh,
                mesh.faceAreas(),
                mesh.cellCentres()
            )
        );

        //- Face based angle
        const scalarField nonOrthoAngle
        (
            radToDeg
            (
149
                Foam::acos(min(scalar(1), faceOrthogonality))
150 151 152 153 154 155 156 157 158 159 160 161
            )
        );

        //- Cell field - max of either face
        volScalarField cellNonOrthoAngle
        (
            IOobject
            (
                "nonOrthoAngle",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
162 163
                IOobject::AUTO_WRITE,
                false
164 165
            ),
            mesh,
166
            dimensionedScalar(dimless, Zero),
167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185
            calculatedFvPatchScalarField::typeName
        );
        //- Take max
        maxFaceToCell(nonOrthoAngle, cellNonOrthoAngle);
        Info<< "    Writing non-orthogonality (angle) to "
            << cellNonOrthoAngle.name() << endl;
        cellNonOrthoAngle.write();
    }

    if (selectedFields.found("faceWeight"))
    {
        volScalarField cellWeights
        (
            IOobject
            (
                "faceWeight",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
186 187
                IOobject::AUTO_WRITE,
                false
188 189
            ),
            mesh,
190
            dimensionedScalar(dimless, Zero),
191 192 193 194 195 196
            wordList                                // wanted bc types
            (
                mesh.boundary().size(),
                calculatedFvPatchScalarField::typeName
            ),
            mesh.weights().boundaryField().types()  // current bc types
197 198
        );
        //- Take min
199
        minFaceToCell(mesh.weights(), cellWeights, false);
200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232
        Info<< "    Writing face interpolation weights (0..0.5) to "
            << cellWeights.name() << endl;
        cellWeights.write();
    }


    // Skewness
    // ~~~~~~~~

    if (selectedFields.found("skewness"))
    {
        //- Face based skewness
        const scalarField faceSkewness
        (
            polyMeshTools::faceSkewness
            (
                mesh,
                mesh.points(),
                mesh.faceCentres(),
                mesh.faceAreas(),
                mesh.cellCentres()
            )
        );

        //- Cell field - max of either face
        volScalarField cellSkewness
        (
            IOobject
            (
                "skewness",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
233 234
                IOobject::AUTO_WRITE,
                false
235 236
            ),
            mesh,
237
            dimensionedScalar(dimless, Zero),
238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263
            calculatedFvPatchScalarField::typeName
        );
        //- Take max
        maxFaceToCell(faceSkewness, cellSkewness);
        Info<< "    Writing face skewness to " << cellSkewness.name() << endl;
        cellSkewness.write();
    }


    // cellDeterminant
    // ~~~~~~~~~~~~~~~

    if (selectedFields.found("cellDeterminant"))
    {
        volScalarField cellDeterminant
        (
            IOobject
            (
                "cellDeterminant",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE,
                false
            ),
            mesh,
264
            dimensionedScalar(dimless, Zero),
265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297
            zeroGradientFvPatchScalarField::typeName
        );
        cellDeterminant.primitiveFieldRef() =
            primitiveMeshTools::cellDeterminant
            (
                mesh,
                mesh.geometricD(),
                mesh.faceAreas(),
                syncTools::getInternalOrCoupledFaces(mesh)
            );
        cellDeterminant.correctBoundaryConditions();
        Info<< "    Writing cell determinant to "
            << cellDeterminant.name() << endl;
        cellDeterminant.write();
    }


    // Aspect ratio
    // ~~~~~~~~~~~~
    if (selectedFields.found("aspectRatio"))
    {
        volScalarField aspectRatio
        (
            IOobject
            (
                "aspectRatio",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE,
                false
            ),
            mesh,
298
            dimensionedScalar(dimless, Zero),
299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321
            zeroGradientFvPatchScalarField::typeName
        );


        scalarField cellOpenness;
        polyMeshTools::cellClosedness
        (
            mesh,
            mesh.geometricD(),
            mesh.faceAreas(),
            mesh.cellVolumes(),
            cellOpenness,
            aspectRatio.ref()
        );

        aspectRatio.correctBoundaryConditions();
        Info<< "    Writing aspect ratio to " << aspectRatio.name() << endl;
        aspectRatio.write();
    }


    // cell type
    // ~~~~~~~~~
322

323 324 325 326 327 328 329 330 331 332 333 334 335 336
    if (selectedFields.found("cellShapes"))
    {
        volScalarField shape
        (
            IOobject
            (
                "cellShapes",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE,
                false
            ),
            mesh,
337
            dimensionedScalar(dimless, Zero),
338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365
            zeroGradientFvPatchScalarField::typeName
        );
        const cellShapeList& cellShapes = mesh.cellShapes();
        forAll(cellShapes, cellI)
        {
            const cellModel& model = cellShapes[cellI].model();
            shape[cellI] = model.index();
        }
        shape.correctBoundaryConditions();
        Info<< "    Writing cell shape (hex, tet etc.) to " << shape.name()
            << endl;
        shape.write();
    }

    if (selectedFields.found("cellVolume"))
    {
        volScalarField V
        (
            IOobject
            (
                "cellVolume",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE,
                false
            ),
            mesh,
366
            dimensionedScalar(dimVolume, Zero),
367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392
            calculatedFvPatchScalarField::typeName
        );
        V.ref() = mesh.V();
        Info<< "    Writing cell volume to " << V.name() << endl;
        V.write();
    }

    if (selectedFields.found("cellVolumeRatio"))
    {
        const scalarField faceVolumeRatio
        (
            polyMeshTools::volRatio
            (
                mesh,
                mesh.V()
            )
        );

        volScalarField cellVolumeRatio
        (
            IOobject
            (
                "cellVolumeRatio",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
393 394
                IOobject::AUTO_WRITE,
                false
395 396
            ),
            mesh,
397
            dimensionedScalar(dimless, Zero),
398 399 400 401 402 403 404 405 406
            calculatedFvPatchScalarField::typeName
        );
        //- Take min
        minFaceToCell(faceVolumeRatio, cellVolumeRatio);
        Info<< "    Writing cell volume ratio to "
            << cellVolumeRatio.name() << endl;
        cellVolumeRatio.write();
    }

407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472
    // minTetVolume
    if (selectedFields.found("minTetVolume"))
    {
        volScalarField minTetVolume
        (
            IOobject
            (
                "minTetVolume",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE,
                false
            ),
            mesh,
            dimensionedScalar("minTetVolume", dimless, GREAT),
            zeroGradientFvPatchScalarField::typeName
        );


        const labelList& own = mesh.faceOwner();
        const labelList& nei = mesh.faceNeighbour();
        const pointField& p = mesh.points();
        forAll(own, facei)
        {
            const face& f = mesh.faces()[facei];
            const point& fc = mesh.faceCentres()[facei];

            {
                const point& ownCc = mesh.cellCentres()[own[facei]];
                scalar& ownVol = minTetVolume[own[facei]];
                forAll(f, fp)
                {
                    scalar tetQual = tetPointRef
                    (
                        p[f[fp]],
                        p[f.nextLabel(fp)],
                        ownCc,
                        fc
                    ).quality();
                    ownVol = min(ownVol, tetQual);
                }
            }
            if (mesh.isInternalFace(facei))
            {
                const point& neiCc = mesh.cellCentres()[nei[facei]];
                scalar& neiVol = minTetVolume[nei[facei]];
                forAll(f, fp)
                {
                    scalar tetQual = tetPointRef
                    (
                        p[f[fp]],
                        p[f.nextLabel(fp)],
                        fc,
                        neiCc
                    ).quality();
                    neiVol = min(neiVol, tetQual);
                }
            }
        }

        minTetVolume.correctBoundaryConditions();
        Info<< "    Writing minTetVolume to " << minTetVolume.name() << endl;
        minTetVolume.write();
    }

473 474 475 476 477 478 479 480 481 482
    if (selectedFields.found("cellRegion"))
    {
        volScalarField cellRegion
        (
            IOobject
            (
                "cellRegion",
                mesh.time().timeName(),
                mesh,
                IOobject::NO_READ,
483 484
                IOobject::AUTO_WRITE,
                false
485 486
            ),
            mesh,
487
            dimensionedScalar(dimless, Zero),
488 489 490 491 492 493 494 495 496 497 498 499
            calculatedFvPatchScalarField::typeName
        );

        regionSplit rs(mesh);
        forAll(rs, celli)
        {
            cellRegion[celli] = rs[celli];
        }
        cellRegion.correctBoundaryConditions();
        Info<< "    Writing cell region to " << cellRegion.name() << endl;
        cellRegion.write();
    }
500

501 502
    Info<< endl;
}