/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2016-2022 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . Application viewFactorsGen Group grpPreProcessingUtilities Description View factors are calculated based on a face agglomeration array (finalAgglom generated by faceAgglomerate utility). Each view factor between the agglomerated faces i and j (Fij) is calculated using a double integral of the sub-areas composing the agglomeration. The patches involved in the view factor calculation are taken from the boundary file and should be part on the group viewFactorWall. ie.: floor { type wall; inGroups 2(wall viewFactorWall); nFaces 100; startFace 3100; } \*---------------------------------------------------------------------------*/ #include "argList.H" #include "fvMesh.H" #include "volFields.H" #include "surfaceFields.H" #include "distributedTriSurfaceMesh.H" #include "meshTools.H" #include "indirectPrimitivePatch.H" #include "DynamicField.H" #include "unitConversion.H" #include "scalarMatrices.H" #include "labelListIOList.H" #include "scalarListIOList.H" #include "singleCellFvMesh.H" #include "IOmapDistribute.H" using namespace Foam; triSurface triangulate ( const polyBoundaryMesh& bMesh, const labelHashSet& includePatches, const labelListIOList& finalAgglom, labelList& triSurfaceToAgglom, const globalIndex& globalNumbering, const polyBoundaryMesh& coarsePatches ) { const polyMesh& mesh = bMesh.mesh(); // Storage for surfaceMesh. Size estimate. DynamicList triangles(mesh.nBoundaryFaces()); label newPatchI = 0; label localTriFaceI = 0; for (const label patchI : includePatches) { const polyPatch& patch = bMesh[patchI]; const pointField& points = patch.points(); label nTriTotal = 0; forAll(patch, patchFaceI) { const face& f = patch[patchFaceI]; faceList triFaces(f.nTriangles(points)); label nTri = 0; f.triangles(points, nTri, triFaces); forAll(triFaces, triFaceI) { const face& f = triFaces[triFaceI]; triangles.append(labelledTri(f[0], f[1], f[2], newPatchI)); nTriTotal++; triSurfaceToAgglom[localTriFaceI++] = globalNumbering.toGlobal ( Pstream::myProcNo(), finalAgglom[patchI][patchFaceI] + coarsePatches[patchI].start() ); } } newPatchI++; } //striSurfaceToAgglom.resize(localTriFaceI-1); triangles.shrink(); // Create globally numbered tri surface triSurface rawSurface(triangles, mesh.points()); // Create locally numbered tri surface triSurface surface ( rawSurface.localFaces(), rawSurface.localPoints() ); // Add patch names to surface surface.patches().setSize(newPatchI); newPatchI = 0; for (const label patchI : includePatches) { const polyPatch& patch = bMesh[patchI]; surface.patches()[newPatchI].index() = patchI; surface.patches()[newPatchI].name() = patch.name(); surface.patches()[newPatchI].geometricType() = patch.type(); newPatchI++; } return surface; } void writeRays ( const fileName& fName, const pointField& compactCf, const pointField& myFc, const labelListList& visibleFaceFaces ) { OFstream str(fName); label vertI = 0; Pout<< "Dumping rays to " << str.name() << endl; forAll(myFc, faceI) { const labelList visFaces = visibleFaceFaces[faceI]; forAll(visFaces, faceRemote) { label compactI = visFaces[faceRemote]; const point& remoteFc = compactCf[compactI]; meshTools::writeOBJ(str, myFc[faceI]); vertI++; meshTools::writeOBJ(str, remoteFc); vertI++; str << "l " << vertI-1 << ' ' << vertI << nl; } } str.flush(); } scalar calculateViewFactorFij ( const vector& i, const vector& j, const vector& dAi, const vector& dAj ) { vector r = i - j; scalar rMag = mag(r); if (rMag > SMALL) { scalar dAiMag = mag(dAi); scalar dAjMag = mag(dAj); vector ni = dAi/dAiMag; vector nj = dAj/dAjMag; scalar cosThetaJ = mag(nj & r)/rMag; scalar cosThetaI = mag(ni & r)/rMag; return ( (cosThetaI*cosThetaJ*dAjMag*dAiMag) /(sqr(rMag)*constant::mathematical::pi) ); } else { return 0; } } void insertMatrixElements ( const globalIndex& globalNumbering, const label fromProcI, const labelListList& globalFaceFaces, const scalarListList& viewFactors, scalarSquareMatrix& matrix ) { forAll(viewFactors, faceI) { const scalarList& vf = viewFactors[faceI]; const labelList& globalFaces = globalFaceFaces[faceI]; label globalI = globalNumbering.toGlobal(fromProcI, faceI); forAll(globalFaces, i) { matrix[globalI][globalFaces[i]] = vf[i]; } } } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { argList::addNote ( "Calculate view factors from face agglomeration array." " The finalAgglom generated by faceAgglomerate utility." ); #include "addRegionOption.H" #include "setRootCase.H" #include "createTime.H" #include "createNamedMesh.H" // Read view factor dictionary IOdictionary viewFactorDict ( IOobject ( "viewFactorsDict", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE ) ); const word viewFactorWall("viewFactorWall"); const bool writeViewFactors = viewFactorDict.getOrDefault("writeViewFactorMatrix", false); const bool dumpRays = viewFactorDict.getOrDefault("dumpRays", false); const label debug = viewFactorDict.getOrDefault