mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: mapFields: add -subtract option to subtract whilst mapping
This commit is contained in:
@ -35,12 +35,13 @@ License
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
template<class Type>
|
||||
template<class Type, class CombineOp>
|
||||
void MapConsistentVolFields
|
||||
(
|
||||
const IOobjectList& objects,
|
||||
const meshToMesh& meshToMeshInterp,
|
||||
const meshToMesh::order& mapOrder
|
||||
const meshToMesh::order& mapOrder,
|
||||
const CombineOp& cop
|
||||
)
|
||||
{
|
||||
const fvMesh& meshSource = meshToMeshInterp.fromMesh();
|
||||
@ -84,7 +85,13 @@ void MapConsistentVolFields
|
||||
);
|
||||
|
||||
// Interpolate field
|
||||
meshToMeshInterp.interpolate(fieldTarget, fieldSource, mapOrder);
|
||||
meshToMeshInterp.interpolate//<Type, eqOp<Type> >
|
||||
(
|
||||
fieldTarget,
|
||||
fieldSource,
|
||||
mapOrder,
|
||||
cop
|
||||
);
|
||||
|
||||
// Write field
|
||||
fieldTarget.write();
|
||||
@ -97,7 +104,12 @@ void MapConsistentVolFields
|
||||
GeometricField<Type, fvPatchField, volMesh> fieldTarget
|
||||
(
|
||||
fieldTargetIOobject,
|
||||
meshToMeshInterp.interpolate(fieldSource, mapOrder)
|
||||
meshToMeshInterp.interpolate//<Type, eqOp<Type> >
|
||||
(
|
||||
fieldSource,
|
||||
mapOrder,
|
||||
cop
|
||||
)
|
||||
);
|
||||
|
||||
// Write field
|
||||
|
||||
263
applications/utilities/preProcessing/mapFields/MapMeshes.H
Normal file
263
applications/utilities/preProcessing/mapFields/MapMeshes.H
Normal file
@ -0,0 +1,263 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
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 <http://www.gnu.org/licenses/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef MapMeshes_H
|
||||
#define MapMeshes_H
|
||||
|
||||
#include "MapVolFields.H"
|
||||
#include "MapConsistentVolFields.H"
|
||||
#include "mapLagrangian.H"
|
||||
#include "UnMapped.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
template<template<class> class CombineOp>
|
||||
void MapConsistentMesh
|
||||
(
|
||||
const fvMesh& meshSource,
|
||||
const fvMesh& meshTarget,
|
||||
const meshToMesh::order& mapOrder
|
||||
)
|
||||
{
|
||||
// Create the interpolation scheme
|
||||
meshToMesh meshToMeshInterp(meshSource, meshTarget);
|
||||
|
||||
Info<< nl
|
||||
<< "Consistently creating and mapping fields for time "
|
||||
<< meshSource.time().timeName() << nl << endl;
|
||||
|
||||
{
|
||||
// Search for list of objects for this time
|
||||
IOobjectList objects(meshSource, meshSource.time().timeName());
|
||||
|
||||
// Map volFields
|
||||
// ~~~~~~~~~~~~~
|
||||
MapConsistentVolFields<scalar>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<scalar>()
|
||||
);
|
||||
MapConsistentVolFields<vector>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<vector>()
|
||||
);
|
||||
MapConsistentVolFields<sphericalTensor>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<sphericalTensor>()
|
||||
);
|
||||
MapConsistentVolFields<symmTensor>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<symmTensor>()
|
||||
);
|
||||
MapConsistentVolFields<tensor>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<tensor>()
|
||||
);
|
||||
}
|
||||
|
||||
{
|
||||
// Search for list of target objects for this time
|
||||
IOobjectList objects(meshTarget, meshTarget.time().timeName());
|
||||
|
||||
// Mark surfaceFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<surfaceScalarField>(objects);
|
||||
UnMapped<surfaceVectorField>(objects);
|
||||
UnMapped<surfaceSphericalTensorField>(objects);
|
||||
UnMapped<surfaceSymmTensorField>(objects);
|
||||
UnMapped<surfaceTensorField>(objects);
|
||||
|
||||
// Mark pointFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<pointScalarField>(objects);
|
||||
UnMapped<pointVectorField>(objects);
|
||||
UnMapped<pointSphericalTensorField>(objects);
|
||||
UnMapped<pointSymmTensorField>(objects);
|
||||
UnMapped<pointTensorField>(objects);
|
||||
}
|
||||
|
||||
mapLagrangian(meshToMeshInterp);
|
||||
}
|
||||
|
||||
|
||||
template<template<class> class CombineOp>
|
||||
void MapSubMesh
|
||||
(
|
||||
const fvMesh& meshSource,
|
||||
const fvMesh& meshTarget,
|
||||
const HashTable<word>& patchMap,
|
||||
const wordList& cuttingPatches,
|
||||
const meshToMesh::order& mapOrder
|
||||
)
|
||||
{
|
||||
// Create the interpolation scheme
|
||||
meshToMesh meshToMeshInterp
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatches
|
||||
);
|
||||
|
||||
Info<< nl
|
||||
<< "Mapping fields for time " << meshSource.time().timeName()
|
||||
<< nl << endl;
|
||||
|
||||
{
|
||||
// Search for list of source objects for this time
|
||||
IOobjectList objects(meshSource, meshSource.time().timeName());
|
||||
|
||||
// Map volFields
|
||||
// ~~~~~~~~~~~~~
|
||||
MapVolFields<scalar>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<scalar>()
|
||||
);
|
||||
MapVolFields<vector>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<vector>()
|
||||
);
|
||||
MapVolFields<sphericalTensor>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<sphericalTensor>()
|
||||
);
|
||||
MapVolFields<symmTensor>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<symmTensor>()
|
||||
);
|
||||
MapVolFields<tensor>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
mapOrder,
|
||||
CombineOp<tensor>()
|
||||
);
|
||||
}
|
||||
|
||||
{
|
||||
// Search for list of target objects for this time
|
||||
IOobjectList objects(meshTarget, meshTarget.time().timeName());
|
||||
|
||||
// Mark surfaceFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<surfaceScalarField>(objects);
|
||||
UnMapped<surfaceVectorField>(objects);
|
||||
UnMapped<surfaceSphericalTensorField>(objects);
|
||||
UnMapped<surfaceSymmTensorField>(objects);
|
||||
UnMapped<surfaceTensorField>(objects);
|
||||
|
||||
// Mark pointFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<pointScalarField>(objects);
|
||||
UnMapped<pointVectorField>(objects);
|
||||
UnMapped<pointSphericalTensorField>(objects);
|
||||
UnMapped<pointSymmTensorField>(objects);
|
||||
UnMapped<pointTensorField>(objects);
|
||||
}
|
||||
|
||||
mapLagrangian(meshToMeshInterp);
|
||||
}
|
||||
|
||||
|
||||
template<template<class> class CombineOp>
|
||||
void MapConsistentSubMesh
|
||||
(
|
||||
const fvMesh& meshSource,
|
||||
const fvMesh& meshTarget,
|
||||
const meshToMesh::order& mapOrder
|
||||
)
|
||||
{
|
||||
HashTable<word> patchMap;
|
||||
HashTable<label> cuttingPatchTable;
|
||||
|
||||
forAll(meshTarget.boundary(), patchi)
|
||||
{
|
||||
if (!isA<processorFvPatch>(meshTarget.boundary()[patchi]))
|
||||
{
|
||||
patchMap.insert
|
||||
(
|
||||
meshTarget.boundary()[patchi].name(),
|
||||
meshTarget.boundary()[patchi].name()
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
cuttingPatchTable.insert
|
||||
(
|
||||
meshTarget.boundaryMesh()[patchi].name(),
|
||||
-1
|
||||
);
|
||||
}
|
||||
}
|
||||
|
||||
MapSubMesh<CombineOp>
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatchTable.toc(),
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -35,12 +35,13 @@ License
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
template<class Type>
|
||||
template<class Type, class CombineOp>
|
||||
void MapVolFields
|
||||
(
|
||||
const IOobjectList& objects,
|
||||
const meshToMesh& meshToMeshInterp,
|
||||
const meshToMesh::order& mapOrder
|
||||
const meshToMesh::order& mapOrder,
|
||||
const CombineOp& cop
|
||||
)
|
||||
{
|
||||
const fvMesh& meshSource = meshToMeshInterp.fromMesh();
|
||||
@ -84,7 +85,13 @@ void MapVolFields
|
||||
);
|
||||
|
||||
// Interpolate field
|
||||
meshToMeshInterp.interpolate(fieldTarget, fieldSource, mapOrder);
|
||||
meshToMeshInterp.interpolate//<Type, eqOp<Type> >
|
||||
(
|
||||
fieldTarget,
|
||||
fieldSource,
|
||||
mapOrder,
|
||||
cop
|
||||
);
|
||||
|
||||
// Write field
|
||||
fieldTarget.write();
|
||||
|
||||
@ -34,11 +34,8 @@ Description
|
||||
|
||||
#include "fvCFD.H"
|
||||
#include "meshToMesh.H"
|
||||
#include "MapVolFields.H"
|
||||
#include "MapConsistentVolFields.H"
|
||||
#include "UnMapped.H"
|
||||
#include "processorFvPatch.H"
|
||||
#include "mapLagrangian.H"
|
||||
#include "MapMeshes.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -46,56 +43,28 @@ void mapConsistentMesh
|
||||
(
|
||||
const fvMesh& meshSource,
|
||||
const fvMesh& meshTarget,
|
||||
const meshToMesh::order& mapOrder
|
||||
const meshToMesh::order& mapOrder,
|
||||
const bool subtract
|
||||
)
|
||||
{
|
||||
// Create the interpolation scheme
|
||||
meshToMesh meshToMeshInterp(meshSource, meshTarget);
|
||||
|
||||
Info<< nl
|
||||
<< "Consistently creating and mapping fields for time "
|
||||
<< meshSource.time().timeName() << nl << endl;
|
||||
|
||||
if (subtract)
|
||||
{
|
||||
// Search for list of objects for this time
|
||||
IOobjectList objects(meshSource, meshSource.time().timeName());
|
||||
|
||||
// Map volFields
|
||||
// ~~~~~~~~~~~~~
|
||||
MapConsistentVolFields<scalar>(objects, meshToMeshInterp, mapOrder);
|
||||
MapConsistentVolFields<vector>(objects, meshToMeshInterp, mapOrder);
|
||||
MapConsistentVolFields<sphericalTensor>
|
||||
MapConsistentMesh<minusEqOp>
|
||||
(
|
||||
objects,
|
||||
meshToMeshInterp,
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder
|
||||
);
|
||||
MapConsistentVolFields<symmTensor>(objects, meshToMeshInterp, mapOrder);
|
||||
MapConsistentVolFields<tensor>(objects, meshToMeshInterp, mapOrder);
|
||||
}
|
||||
|
||||
else
|
||||
{
|
||||
// Search for list of target objects for this time
|
||||
IOobjectList objects(meshTarget, meshTarget.time().timeName());
|
||||
|
||||
// Mark surfaceFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<surfaceScalarField>(objects);
|
||||
UnMapped<surfaceVectorField>(objects);
|
||||
UnMapped<surfaceSphericalTensorField>(objects);
|
||||
UnMapped<surfaceSymmTensorField>(objects);
|
||||
UnMapped<surfaceTensorField>(objects);
|
||||
|
||||
// Mark pointFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<pointScalarField>(objects);
|
||||
UnMapped<pointVectorField>(objects);
|
||||
UnMapped<pointSphericalTensorField>(objects);
|
||||
UnMapped<pointSymmTensorField>(objects);
|
||||
UnMapped<pointTensorField>(objects);
|
||||
MapConsistentMesh<eqOp>
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
|
||||
mapLagrangian(meshToMeshInterp);
|
||||
}
|
||||
|
||||
|
||||
@ -105,57 +74,32 @@ void mapSubMesh
|
||||
const fvMesh& meshTarget,
|
||||
const HashTable<word>& patchMap,
|
||||
const wordList& cuttingPatches,
|
||||
const meshToMesh::order& mapOrder
|
||||
const meshToMesh::order& mapOrder,
|
||||
const bool subtract
|
||||
)
|
||||
{
|
||||
// Create the interpolation scheme
|
||||
meshToMesh meshToMeshInterp
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatches
|
||||
);
|
||||
|
||||
Info<< nl
|
||||
<< "Mapping fields for time " << meshSource.time().timeName()
|
||||
<< nl << endl;
|
||||
|
||||
if (subtract)
|
||||
{
|
||||
// Search for list of source objects for this time
|
||||
IOobjectList objects(meshSource, meshSource.time().timeName());
|
||||
|
||||
// Map volFields
|
||||
// ~~~~~~~~~~~~~
|
||||
MapVolFields<scalar>(objects, meshToMeshInterp, mapOrder);
|
||||
MapVolFields<vector>(objects, meshToMeshInterp, mapOrder);
|
||||
MapVolFields<sphericalTensor>(objects, meshToMeshInterp, mapOrder);
|
||||
MapVolFields<symmTensor>(objects, meshToMeshInterp, mapOrder);
|
||||
MapVolFields<tensor>(objects, meshToMeshInterp, mapOrder);
|
||||
MapSubMesh<minusEqOp>
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatches,
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
|
||||
else
|
||||
{
|
||||
// Search for list of target objects for this time
|
||||
IOobjectList objects(meshTarget, meshTarget.time().timeName());
|
||||
|
||||
// Mark surfaceFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<surfaceScalarField>(objects);
|
||||
UnMapped<surfaceVectorField>(objects);
|
||||
UnMapped<surfaceSphericalTensorField>(objects);
|
||||
UnMapped<surfaceSymmTensorField>(objects);
|
||||
UnMapped<surfaceTensorField>(objects);
|
||||
|
||||
// Mark pointFields as unmapped
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
UnMapped<pointScalarField>(objects);
|
||||
UnMapped<pointVectorField>(objects);
|
||||
UnMapped<pointSphericalTensorField>(objects);
|
||||
UnMapped<pointSymmTensorField>(objects);
|
||||
UnMapped<pointTensorField>(objects);
|
||||
MapSubMesh<eqOp>
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatches,
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
|
||||
mapLagrangian(meshToMeshInterp);
|
||||
}
|
||||
|
||||
|
||||
@ -163,40 +107,28 @@ void mapConsistentSubMesh
|
||||
(
|
||||
const fvMesh& meshSource,
|
||||
const fvMesh& meshTarget,
|
||||
const meshToMesh::order& mapOrder
|
||||
const meshToMesh::order& mapOrder,
|
||||
const bool subtract
|
||||
)
|
||||
{
|
||||
HashTable<word> patchMap;
|
||||
HashTable<label> cuttingPatchTable;
|
||||
|
||||
forAll(meshTarget.boundary(), patchi)
|
||||
if (subtract)
|
||||
{
|
||||
if (!isA<processorFvPatch>(meshTarget.boundary()[patchi]))
|
||||
{
|
||||
patchMap.insert
|
||||
(
|
||||
meshTarget.boundary()[patchi].name(),
|
||||
meshTarget.boundary()[patchi].name()
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
cuttingPatchTable.insert
|
||||
(
|
||||
meshTarget.boundaryMesh()[patchi].name(),
|
||||
-1
|
||||
);
|
||||
}
|
||||
MapConsistentSubMesh<minusEqOp>
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
MapConsistentSubMesh<eqOp>
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
|
||||
mapSubMesh
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatchTable.toc(),
|
||||
mapOrder
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
@ -288,6 +220,11 @@ int main(int argc, char *argv[])
|
||||
"word",
|
||||
"specify the mapping method"
|
||||
);
|
||||
argList::addBoolOption
|
||||
(
|
||||
"subtract",
|
||||
"subtract mapped source from target"
|
||||
);
|
||||
|
||||
argList args(argc, argv);
|
||||
|
||||
@ -350,6 +287,13 @@ int main(int argc, char *argv[])
|
||||
Info<< "Mapping method: " << mapMethod << endl;
|
||||
}
|
||||
|
||||
const bool subtract = args.optionFound("subtract");
|
||||
if (subtract)
|
||||
{
|
||||
Info<< "Subtracting mapped source field from target" << endl;
|
||||
}
|
||||
|
||||
|
||||
#include "createTimes.H"
|
||||
|
||||
HashTable<word> patchMap;
|
||||
@ -431,7 +375,13 @@ int main(int argc, char *argv[])
|
||||
|
||||
if (consistent)
|
||||
{
|
||||
mapConsistentSubMesh(meshSource, meshTarget, mapOrder);
|
||||
mapConsistentSubMesh
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
@ -441,7 +391,8 @@ int main(int argc, char *argv[])
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatches,
|
||||
mapOrder
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
}
|
||||
@ -503,7 +454,13 @@ int main(int argc, char *argv[])
|
||||
|
||||
if (consistent)
|
||||
{
|
||||
mapConsistentSubMesh(meshSource, meshTarget, mapOrder);
|
||||
mapConsistentSubMesh
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
@ -513,7 +470,8 @@ int main(int argc, char *argv[])
|
||||
meshTarget,
|
||||
patchMap,
|
||||
addProcessorPatches(meshTarget, cuttingPatches),
|
||||
mapOrder
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
}
|
||||
@ -629,7 +587,8 @@ int main(int argc, char *argv[])
|
||||
(
|
||||
meshSource,
|
||||
meshTarget,
|
||||
mapOrder
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
else
|
||||
@ -640,7 +599,8 @@ int main(int argc, char *argv[])
|
||||
meshTarget,
|
||||
patchMap,
|
||||
addProcessorPatches(meshTarget, cuttingPatches),
|
||||
mapOrder
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
}
|
||||
@ -679,7 +639,7 @@ int main(int argc, char *argv[])
|
||||
|
||||
if (consistent)
|
||||
{
|
||||
mapConsistentMesh(meshSource, meshTarget, mapOrder);
|
||||
mapConsistentMesh(meshSource, meshTarget, mapOrder, subtract);
|
||||
}
|
||||
else
|
||||
{
|
||||
@ -689,7 +649,8 @@ int main(int argc, char *argv[])
|
||||
meshTarget,
|
||||
patchMap,
|
||||
cuttingPatches,
|
||||
mapOrder
|
||||
mapOrder,
|
||||
subtract
|
||||
);
|
||||
}
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user