/// <summary> /// Reinit on un-cut cells. /// </summary> /// <param name="Phi">The level set</param> /// <param name="ReInitSpecies">Cell mask which is to be reinitialized</param> /// <param name="sign">Sign of the level set for this <paramref name="ReInitSpecies"/></param> /// <param name="_Accepted">CellMask which is taken as boundary values</param> /// <param name="GradPhi">Level Set gradient</param> /// <param name="callBack">A delegate, which might be called after the execution of the reinitialization</param> public void Reinitialize(SinglePhaseField Phi, CellMask ReInitSpecies, double sign, CellMask _Accepted, //ConventionalDGField[] ExtProperty, double[][] ExtPropertyMin, double[][] ExtPropertyMax, VectorField <SinglePhaseField> GradPhi, Action <int> callBack) // { using (new FuncTrace()) { Tracer.InstrumentationSwitch = false; // lots of tracing on calls acting on singe cells causes massive overhead (up to 5x slower). Stpw_total.Start(); SinglePhaseField DiffusionCoeff = new SinglePhaseField(new Basis(this.GridDat, 1), "DiffusionCoeff"); // check args and init // =================== /* * ExtVelSolver extVelSlv = null; * if(ExtProperty != null) { * if(ExtProperty.Length != ExtPropertyMin.Length) * throw new ArgumentException(); * if(ExtProperty.Length != ExtPropertyMax.Length) * throw new ArgumentException(); * * extVelSlv = new ExtVelSolver(ExtProperty[0].Basis); * } */ BitArray Acceped_Mutuable = _Accepted.GetBitMask().CloneAs(); BitArray Trial_Mutuable = ((_Accepted.AllNeighbourCells().Intersect(ReInitSpecies)).Except(_Accepted)).GetBitMask().CloneAs(); BitArray Recalc_Mutuable = Trial_Mutuable.CloneAs(); BitArray PosSpecies_Bitmask = ReInitSpecies.GetBitMask(); int J = this.GridDat.Cells.Count; int D = this.GridDat.SpatialDimension; int N = this.LevelSetBasis.Length; double _sign = sign >= 0 ? 1.0 : -1.0; double[] PhiAvg = m_PhiAvg; if (PhiAvg == null) { throw new ApplicationException(); } foreach (int jCell in _Accepted.ItemEnum) { PhiAvg[jCell] = Phi.GetMeanValue(jCell); } int NoOfNew; { var Neu = ReInitSpecies.Except(_Accepted); NoOfNew = Neu.NoOfItemsLocally; Phi.Clear(Neu); Phi.AccConstant(_sign, Neu); foreach (int jCell in Neu.ItemEnum) { PhiAvg[jCell] = 1.0e10; } } if (this.GridDat.MpiSize > 1) { throw new NotSupportedException("Currently not MPI parallel."); } for (int d = 0; d < this.GridDat.SpatialDimension; d++) { if (!GradPhi[d].Basis.Equals(Phi.Basis)) { throw new ArgumentException("Level-set and level-set gradient field should have the same DG basis."); // ein grad niedriger wrürde auch genügen... } } // perform marching... // =================== // update gradient for cut-cells GradPhi.Clear(_Accepted); GradPhi.Gradient(1.0, Phi, _Accepted); // marching loop../ int cnt = 0; while (true) { cnt++; CellMask Recalc = new CellMask(this.GridDat, Recalc_Mutuable); CellMask Accepted = new CellMask(this.GridDat, Acceped_Mutuable); CellMask Trial = new CellMask(this.GridDat, Trial_Mutuable); int NoOfTrial = Trial.NoOfItemsLocally; int NoOfAccpt = Accepted.NoOfItemsLocally; int NoOfRcalc = Recalc.NoOfItemsLocally; if (Trial.NoOfItemsLocally <= 0) { //Ploti(Recalc, Accepted, Trial, Phi, Phi_gradient, optEikonalOut, cnt); break; } // Local solver for all 'Recalc'-cells // -------------------------------------- if (Recalc.NoOfItemsLocally > 0) { this.LocalSolve(Accepted, Recalc, Phi, GradPhi, _sign, DiffusionCoeff); } // find the next cell to accept // ---------------------------- // get mean value in all cells foreach (int jCell in Recalc.ItemEnum) { PhiAvg[jCell] = Phi.GetMeanValue(jCell); Recalc_Mutuable[jCell] = false; } //Ploti(Recalc, Accepted, Trial, Phi, Phi_gradient, optEikonalOut, cnt); // find trial-cell with minimum average value // this should be done with heap-sort (see fast-marching algorithm) int jCellAccpt = int.MaxValue; double TrialMin = double.MaxValue; foreach (int jCell in Trial.ItemEnum) { if (PhiAvg[jCell] * _sign < TrialMin) { TrialMin = PhiAvg[jCell] * _sign; jCellAccpt = jCell; } } if (callBack != null) { callBack(cnt); } /* * // update the gradient * // ------------------- * * this.Stpw_gradientEval.Start(); * gradModule.GradientUpdate(jCellAccpt, Acceped_Mutuable, Phi, GradPhi); * this.Stpw_gradientEval.Stop(); * * /* * // solve for the extension properties * // ---------------------------------- * * if(ExtProperty != null) { * int[] Neight, dummy33; * GridDat.Cells.GetCellNeighbours(jCellAccpt, GridData.CellData.GetCellNeighbours_Mode.ViaEdges, out Neight, out dummy33); * * for(int iComp = 0; iComp < ExtProperty.Length; iComp++) { * * ExtPropertyMax[iComp][jCellAccpt] = -double.MaxValue; * ExtPropertyMin[iComp][jCellAccpt] = double.MaxValue; * * foreach(int jNeig in Neight) { * if(Acceped_Mutuable[jNeig]) { * ExtPropertyMax[iComp][jCellAccpt] = Math.Max(ExtPropertyMax[iComp][jCellAccpt], ExtPropertyMax[iComp][jNeig]); * ExtPropertyMin[iComp][jCellAccpt] = Math.Min(ExtPropertyMin[iComp][jCellAccpt], ExtPropertyMin[iComp][jNeig]); * } * } * * this.Stpw_extVelSolver.Start(); * extVelSlv.ExtVelSolve_Far(Phi, GradPhi, ExtProperty[iComp], ref ExtPropertyMin[iComp][jCellAccpt], ref ExtPropertyMax[iComp][jCellAccpt], jCellAccpt, Accepted, _sign); * this.Stpw_extVelSolver.Stop(); * } * } * /* * { * int[] Neight, dummy33; * GridDat.Cells.GetCellNeighbours(jCellAccpt, GridData.CellData.GetCellNeighbours_Mode.ViaEdges, out Neight, out dummy33); * foreach(int jNeig in Neight) { * if(Acceped_Mutuable[jNeig]) { * plotDependencyArrow(cnt, jCellAccpt, jNeig); * } * } * } */ // the mimium is moved to accepted // ------------------------------- Acceped_Mutuable[jCellAccpt] = true; Trial_Mutuable[jCellAccpt] = false; Recalc_Mutuable[jCellAccpt] = false; NoOfNew--; // recalc on all neighbours // ------------------------ int[] Neighs, dummy; this.GridDat.GetCellNeighbours(jCellAccpt, GetCellNeighbours_Mode.ViaEdges, out Neighs, out dummy); foreach (int jNeig in Neighs) { if (!Acceped_Mutuable[jNeig] && PosSpecies_Bitmask[jNeig]) { Trial_Mutuable[jNeig] = true; Recalc_Mutuable[jNeig] = true; } } } if (NoOfNew > 0) { throw new ArithmeticException("Unable to perform reinitialization for all requested cells - maybe they are not reachable from the initialy 'accepted' domain?"); } //PlottAlot("dependencies.csv"); Tracer.InstrumentationSwitch = true; Stpw_total.Stop(); } }
/// <summary> /// /// </summary> /// <param name="InitialMass"></param> /// <param name="Temperature"></param> /// <returns></returns> public override double GetMassDeterminedThermodynamicPressure(double InitialMass, SinglePhaseField Temperature) { throw new NotImplementedException(); }
void TestAgglomeration_Projection(int quadOrder, MultiphaseCellAgglomerator Agg) { var Bmask = LsTrk.Regions.GetSpeciesMask("B"); var BbitMask = Bmask.GetBitMask(); var XDGmetrics = LsTrk.GetXDGSpaceMetrics(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, 1); int degree = Math.Max(0, this.u.Basis.Degree - 1); _2D SomePolynomial; if (degree == 0) { SomePolynomial = (x, y) => 3.2; } else if (degree == 1) { SomePolynomial = (x, y) => 3.2 + x - 2.4 * y; } else { SomePolynomial = (x, y) => x * x + y * y; } // ------------------------------------------------ // project the polynomial onto a single-phase field // ------------------------------------------------ SinglePhaseField NonAgglom = new SinglePhaseField(new Basis(this.GridData, degree), "NonAgglom"); NonAgglom.ProjectField(1.0, SomePolynomial.Vectorize(), new CellQuadratureScheme(true, Bmask)); // ------------------------------------------------ // project the polynomial onto the aggomerated XDG-field // ------------------------------------------------ var qsh = XDGmetrics.XQuadSchemeHelper; // Compute the inner product of 'SomePolynomial' with the cut-cell-basis, ... SinglePhaseField xt = new SinglePhaseField(NonAgglom.Basis, "test"); xt.ProjectField(1.0, SomePolynomial.Vectorize(), qsh.GetVolumeQuadScheme(this.LsTrk.GetSpeciesId("B")).Compile(this.GridData, Agg.CutCellQuadratureOrder)); CoordinateMapping map = xt.Mapping; // ... change to the cut-cell-basis, ... double[] xVec = xt.CoordinateVector.ToArray(); Agg.ManipulateMatrixAndRHS(default(MsrMatrix), xVec, map, null); { double[] xVec2 = xVec.CloneAs(); Agg.ClearAgglomerated(xVec2, map); // clearing should have no effect now. double Err3 = GenericBlas.L2DistPow2(xVec, xVec2).MPISum().Sqrt(); Assert.LessOrEqual(Err3, 0.0); } // ... multiply with the inverse mass-matrix, ... var Mfact = XDGmetrics.MassMatrixFactory; BlockMsrMatrix MassMtx = Mfact.GetMassMatrix(map, false); Agg.ManipulateMatrixAndRHS(MassMtx, default(double[]), map, map); var InvMassMtx = MassMtx.InvertBlocks(OnlyDiagonal: true, Subblocks: true, ignoreEmptyBlocks: true, SymmetricalInversion: true); double[] xAggVec = new double[xVec.Length]; InvMassMtx.SpMV(1.0, xVec, 0.0, xAggVec); // ... and extrapolate. // Since we projected a polynomial, the projection is exact and after extrapolation, must be equal // to the polynomial. xt.CoordinateVector.SetV(xAggVec); Agg.Extrapolate(xt.Mapping); // ------------------------------------------------ // check the error // ------------------------------------------------ //double Err2 = xt.L2Error(SomePolynomial.Vectorize(), new CellQuadratureScheme(domain: Bmask)); double Err2 = NonAgglom.L2Error(xt); Console.WriteLine("Agglom: projection and extrapolation error: " + Err2); Assert.LessOrEqual(Err2, 1.0e-10); }
/// <summary> /// Based on the Ideas by /// C. Basting and D. Kuzmin, /// “A minimization-based finite element formulation for interface-preserving level set reinitialization”, /// Computing, vol. 95, no. 1, pp. 13–25, Dec. 2012. /// Create Spatial Operators and build the corresponding Matrices /// For the Left-Hand Side of the ReInitProblem /// RHS is computed on the fly in <see cref="ReInitSingleStep"/> /// The Bulk component is constant unless the grid changes, thus it is computed in <see cref="BuildOperators(CellQuadratureScheme)"/>. /// The Interface component changes with its motion. /// This component is calculated in <see cref="UpdateOperators(CellQuadratureScheme)"/>. /// </summary> /// <param name="LSTrck"></param> /// <param name="Control">various parameters <paramref name="EllipticReinitControl"/></param> /// <param name="HMFOrder">order of tghe interface quadrature</param> public EllipticReInit(LevelSetTracker LSTrck, EllipticReInitAlgoControl Control, SinglePhaseField LevelSetForReInit = null) { this.Control = Control; this.LevelSetTracker = LSTrck; if (LevelSetForReInit == null) { Phi = LevelSetTracker.LevelSets[0] as SinglePhaseField; } else { Phi = LevelSetForReInit; } this.underrelaxation = Control.underrelaxation; Residual = new SinglePhaseField(Phi.Basis); OldPhi = new SinglePhaseField(Phi.Basis); NewPhi = new SinglePhaseField(Phi.Basis); foreach (SinglePhaseField f in new List <SinglePhaseField> { Residual, OldPhi, NewPhi }) { f.Clear(); f.Acc(1.0, Phi); } this.D = LevelSetTracker.GridDat.SpatialDimension; this.ConvergenceCriterion = Control.ConvergenceCriterion; this.MaxIteration = Control.MaxIt; double PenaltyBase = ((double)((Phi.Basis.Degree + 1) * (Phi.Basis.Degree + D))) / ((double)D); /// Choose Forms according to Upwinding or Central Fluxes string[] paramNames; int noOfParamFields; IEquationComponent BulkForm; RHSForm myRHSForm; LevelSetGradient = new VectorField <SinglePhaseField>(D, Phi.Basis, "LevelSetGradient", SinglePhaseField.Factory); MeanLevelSetGradient = new VectorField <SinglePhaseField>(D, new Basis(Phi.GridDat, 0), "MeanLevelSetGradient", SinglePhaseField.Factory); if (Control.Upwinding) { paramNames = new string[] { "OldLevelSet", "MeanLevelSetGradient[0]", "MeanLevelSetGradient[1]" }; noOfParamFields = D; LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient); parameterFields = ArrayTools.Cat(new SinglePhaseField[] { OldPhi }, MeanLevelSetGradient.ToArray()); //throw new NotImplementedException("ToDO"); BulkForm = new EllipticReInitUpwindForm_Laplace(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck); myRHSForm = new EllipticReInitUpwindForm_RHS(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck); OldDirection = new double[MeanLevelSetGradient.CoordinateVector.ToArray().Length]; for (int i = 0; i < MeanLevelSetGradient.CoordinateVector.Length; i++) { OldDirection[i] = Math.Sign(MeanLevelSetGradient.CoordinateVector[i]); } NewDirection = OldDirection.CloneAs(); } else { paramNames = new string[] { }; noOfParamFields = 0; parameterFields = new SinglePhaseField[] { }; BulkForm = new CentralDifferencesLHSForm(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck.GridDat.Cells.cj); myRHSForm = new CentralDifferencesRHSForm(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck); } // SIP for the bulk Phase //this.Operator_bulk = new SpatialOperator(1, noOfParamFields, 1, QuadOrderFunc.SumOfMaxDegrees(1, RoundUp: false), variableNames); this.Operator_bulk = BulkForm.Operator(); // Zero at the Interface // Calculate Quadrature Order Func <int[], int[], int[], int> InterfaceQuadOrder; InterfaceQuadOrder = QuadOrderFunc.FixedOrder(Phi.Basis.Degree * 2 + 2); // Generate Interface Operator this.Operator_interface = (new EllipticReInitInterfaceForm(Control.PenaltyMultiplierInterface * PenaltyBase, LSTrck)).XOperator(InterfaceQuadOrder); // Nonlinear Part on the RHS // switch for the potential functions switch (Control.Potential) { case ReInitPotential.BastingDoubleWell: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.DoubleWell(d, b)); break; }; case ReInitPotential.BastingSingleWell: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.SingleWell(d, b)); break; }; case ReInitPotential.SingleWellNear: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.SingleWellNear(d, b)); break; }; case ReInitPotential.P4DoubleWell: { Console.WriteLine("Warning - This Option for Elliptic ReInit does not work well"); myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.DoubleWellAlternative(d, b)); break; }; case ReInitPotential.SingleWellOnCutDoubleWellElse: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.SingleWellOnCutDoubleWellElse(d, b)); break; } } Operator_RHS = myRHSForm.Operator(QuadOrderFunc.SumOfMaxDegrees(2, RoundUp: true)); // The result of the nonlinear part on the rhs is projected on a single-phase field RHSField = new SinglePhaseField(Phi.Basis, "RHS"); OpMatrix = new MsrMatrix(this.Phi.Mapping, this.Phi.Mapping); OpAffine = new double[OpMatrix.RowPartitioning.LocalLength]; // Matrix and RHS for the Bulk component OpMatrix_bulk = new MsrMatrix(this.Phi.Mapping, this.Phi.Mapping); OpAffine_bulk = new double[OpMatrix.RowPartitioning.LocalLength]; // Matrix and RHS for the Interface Penalty OpMatrix_interface = new MsrMatrix(this.Phi.Mapping, this.Phi.Mapping); OpAffine_interface = new double[OpMatrix.RowPartitioning.LocalLength]; // Init Parameter Fields OldPhi.Clear(); OldPhi.Acc(1.0, Phi); // Compute Matrices UpdateBulkMatrix(); }
/// <summary> /// DG field instantiation /// </summary> protected override void CreateFields() { base.CreateFields(); ResiualKP1 = new SinglePhaseField(new Basis(this.GridData, T.Basis.Degree + 1), "ResidualKP1"); base.IOFields.Add(ResiualKP1); }
public void MakeContinuous(SinglePhaseField DGLevelSet, SinglePhaseField LevelSet, CellMask Domain) { CDGField.ProjectDGField(1.0, DGLevelSet, Domain); LevelSet.Clear(); CDGField.AccToDGField(1.0, LevelSet); }
/// <summary> /// Ctor. /// </summary> /// <param name="solverConf"></param> /// <param name="sparseSolver"></param> /// <param name="MatAsmblyCorrector"></param> /// <param name="VelocityDivergence"></param> /// <param name="Velocity_Intrmed"></param> /// <param name="DivB4"></param> /// <param name="BDF"></param> /// <param name="Temperature"></param> /// <param name="EoS"></param> public LowMachSolverCorrector(SolverConfiguration solverConf, ISparseSolver sparseSolver, SIMPLEMatrixAssembly MatAsmblyCorrector, SIMPLEOperator[] VelocityDivergence, VectorField <SinglePhaseField> Velocity_Intrmed, SinglePhaseField DivB4, BDFScheme BDF, ScalarFieldHistory <SinglePhaseField> Temperature, MaterialLaw EoS, double[] RHSManuDivKontiOperatorAffine = null) : base(solverConf, sparseSolver) { this.MatAsmblyCorrector = MatAsmblyCorrector; this.VelocityDivergence = VelocityDivergence; this.Velocity_Intrmed = Velocity_Intrmed; this.DivB4 = DivB4; this.BDF = BDF; this.Temperature = Temperature; this.EoS = EoS; this.RHSManuDivKontiOperatorAffine = RHSManuDivKontiOperatorAffine; }
/// <summary> /// Ctor without pressure stabilization. /// </summary> /// <param name="solverConfig"></param> /// <param name="sparseSolver"></param> /// <param name="VelocityDivergence"></param> /// <param name="MatAsmblyCorrector"></param> /// <param name="Velocity_Intrmed"></param> /// <param name="DivB4"></param> public MultiphaseSolverCorrector(SolverConfiguration solverConfig, ISparseSolver sparseSolver, SIMPLEOperator[] VelocityDivergence, SIMPLEMatrixAssembly MatAsmblyCorrector, VectorField <SinglePhaseField> Velocity_Intrmed, SinglePhaseField DivB4) : base(solverConfig, sparseSolver) { if ((VelocityDivergence.Length != solverConfig.SpatialDimension) || (Velocity_Intrmed.Dim != solverConfig.SpatialDimension)) { throw new ArgumentException("Mismatch of dimensions!"); } m_VelocityDivergence = VelocityDivergence; m_MatAsmblyCorrector = MatAsmblyCorrector; m_Velocity_Intrmed = Velocity_Intrmed; m_DivB4 = DivB4; }
/// <summary> /// Calculates the lift and drag force on given edges, e.g edges around an airfoil. /// Currently only in 2D!!! /// </summary> /// <param name="edgeTagName">EdgeTag on which the drag force will be calculated</param> /// <returns>total drag force</returns> /// <remarks>It assumes that pressure and velocity fields exists</remarks> public static Query LiftAndDragForce(String edgeTagName) { return(delegate(IApplication <AppControl> app, double time) { if (dragIntegral == null) { densityField = app.IOFields.SingleOrDefault((f) => f.Identification == "rho"); m0Field = app.IOFields.SingleOrDefault((f) => f.Identification == "m0"); m1Field = app.IOFields.SingleOrDefault((f) => f.Identification == "m1"); energyField = app.IOFields.SingleOrDefault((f) => f.Identification == "rhoE"); DrhoDx = new SinglePhaseField(densityField.Basis); DrhoDy = new SinglePhaseField(densityField.Basis); Dm0Dx = new SinglePhaseField(m0Field.Basis); Dm0Dy = new SinglePhaseField(m0Field.Basis); Dm1Dx = new SinglePhaseField(m1Field.Basis); Dm1Dy = new SinglePhaseField(m1Field.Basis); CNSControl c = app.Control as CNSControl; Program prog = app as Program; byte edgeTag = app.Grid.EdgeTagNames.First(item => item.Value.Equals(edgeTagName)).Key; CoordinateMapping mapping = new CoordinateMapping( densityField, m0Field, m1Field, energyField, DrhoDx, DrhoDy, Dm0Dx, Dm0Dy, Dm1Dx, Dm1Dy); dragIntegral = new EdgeIntegral(app.GridData, edgeTag, new ForceFlux( c.ReynoldsNumber, prog.SpeciesMap.GetMaterial(double.NaN), 0), mapping); liftIntegral = new EdgeIntegral(app.GridData, edgeTag, new ForceFlux( c.ReynoldsNumber, prog.SpeciesMap.GetMaterial(double.NaN), 1), mapping); if (logger == null && app.CurrentSessionInfo.ID != Guid.Empty && app.MPIRank == 0) { logger = app.DatabaseDriver.FsDriver.GetNewLog("LiftAndDragForce", app.CurrentSessionInfo.ID); string header = "PhysTime\t LiftForce\t DragForce"; logger.WriteLine(header); } } DrhoDx.Clear(); DrhoDy.Clear(); DrhoDx.Derivative(1.0, densityField, 0); DrhoDy.Derivative(1.0, densityField, 1); Dm0Dx.Clear(); Dm0Dy.Clear(); Dm1Dx.Clear(); Dm1Dy.Clear(); Dm0Dx.Derivative(1.0, m0Field, 0); Dm0Dy.Derivative(1.0, m0Field, 1); Dm1Dx.Derivative(1.0, m1Field, 0); Dm1Dy.Derivative(1.0, m1Field, 1); double lift = liftIntegral.Evaluate(); double drag = dragIntegral.Evaluate(); if (logger != null && app.MPIRank == 0) { string line = time + "\t" + lift + "\t" + drag; logger.WriteLine(line); logger.Flush(); } return drag; }); }
/// <summary> /// Algorithm to find an inflection point along a curve in the direction of the gradient /// </summary> /// <param name="gridData">The corresponding grid</param> /// <param name="field">The DG field, which shall be evalauted</param> /// <param name="results"> /// The points along the curve (first point has to be user defined) /// Lenghts --> [0]: maxIterations + 1, [1]: 5 /// [1]: x | y | function values | second derivatives | step sizes /// </param> /// <param name="iterations">The amount of iterations needed in order to find the inflection point</param> /// <param name="converged">Has an inflection point been found?</param> /// <param name = "jLocal" >Local cell index of the inflection point</ param > /// <param name="byFlux">Option for the calculation of the first and second order derivatives, default: true</param> private void WalkOnCurve(GridData gridData, SinglePhaseField field, MultidimensionalArray results, out int iterations, out bool converged, out int jLocal, bool byFlux = true) { // Init converged = false; // Current (global) point double[] currentPoint = new double[] { results[0, 0], results[0, 1] }; // Compute global cell index of current point gridData.LocatePoint(currentPoint, out long GlobalId, out long GlobalIndex, out bool IsInside, out bool OnThisProcess); // Compute local node set NodeSet nodeSet = ShockFindingExtensions.GetLocalNodeSet(gridData, currentPoint, (int)GlobalIndex); // Get local cell index of current point int j0Grd = gridData.CellPartitioning.i0; jLocal = (int)(GlobalIndex - j0Grd); // Evaluate the second derivative try { if (byFlux) { results[0, 3] = SecondDerivativeByFlux(field, jLocal, nodeSet); } else { results[0, 3] = ShockFindingExtensions.SecondDerivative(field, jLocal, nodeSet); } } catch (NotSupportedException) { iterations = 0; return; } // Evaluate the function MultidimensionalArray f = MultidimensionalArray.Create(1, 1); field.Evaluate(jLocal, 1, nodeSet, f); results[0, 2] = f[0, 0]; // Set initial step size to 0.5 * h_minGlobal results[0, 4] = 0.5 * gridData.Cells.h_minGlobal; int n = 1; while (n < results.Lengths[0]) { // Evaluate the gradient of the current point double[] gradient; if (byFlux) { gradient = NormalizedGradientByFlux(field, jLocal, nodeSet); } else { gradient = ShockFindingExtensions.NormalizedGradient(field, jLocal, nodeSet); } // Compute new point along curve currentPoint[0] = currentPoint[0] + gradient[0] * results[n - 1, 4]; currentPoint[1] = currentPoint[1] + gradient[1] * results[n - 1, 4]; // New point has been calculated --> old node set is invalid nodeSet = null; // Check if new point is still in the same cell or has moved to one of its neighbours if (!gridData.Cells.IsInCell(currentPoint, jLocal)) { // Get indices of cell neighbours gridData.GetCellNeighbours(jLocal, GetCellNeighbours_Mode.ViaVertices, out int[] cellNeighbours, out int[] connectingEntities); double[] newLocalCoord = new double[currentPoint.Length]; bool found = false; // Find neighbour foreach (int neighbour in cellNeighbours) { if (gridData.Cells.IsInCell(currentPoint, neighbour, newLocalCoord)) { // If neighbour has been found, update jLocal = neighbour + j0Grd; nodeSet = ShockFindingExtensions.GetLocalNodeSet(gridData, currentPoint, jLocal); found = true; break; } } if (found == false) { iterations = n; return; } } else { // New point is still in the same cell --> update only the coordiantes (local node set) nodeSet = ShockFindingExtensions.GetLocalNodeSet(gridData, currentPoint, jLocal + j0Grd); } // Update output results[n, 0] = currentPoint[0]; results[n, 1] = currentPoint[1]; // Evaluate the function f.Clear(); field.Evaluate(jLocal, 1, nodeSet, f); results[n, 2] = f[0, 0]; // Evaluate the second derivative of the new point try { if (byFlux) { results[n, 3] = SecondDerivativeByFlux(field, jLocal, nodeSet); } else { results[n, 3] = ShockFindingExtensions.SecondDerivative(field, jLocal, nodeSet); } } catch (NotSupportedException) { break; } // Check if sign of second derivative has changed // Halve the step size and change direction if (Math.Sign(results[n, 3]) != Math.Sign(results[n - 1, 3])) { results[n, 4] = -results[n - 1, 4] / 2; } else { results[n, 4] = results[n - 1, 4]; } n++; // Termination criterion // Remark: n has already been incremented! double[] diff = new double[currentPoint.Length]; diff[0] = results[n - 1, 0] - results[n - 2, 0]; diff[1] = results[n - 1, 1] - results[n - 2, 1]; if (diff.L2Norm() < 1e-12) { converged = true; break; } } iterations = n; return; }
void DerivativeByFluxLinear(SinglePhaseField fin, SinglePhaseField fres, int d, SinglePhaseField fBnd) { var Op = (new LinearDerivFlx(d)).Operator(); MsrMatrix OpMtx = new MsrMatrix(fres.Mapping, fin.Mapping); double[] OpAff = new double[fres.Mapping.LocalLength]; Op.ComputeMatrixEx(fin.Mapping, new DGField[] { fBnd }, fres.Mapping, OpMtx, OpAff, OnlyAffine: false); fres.Clear(); fres.CoordinateVector.Acc(1.0, OpAff); OpMtx.SpMVpara(1.0, fin.CoordinateVector, 1.0, fres.CoordinateVector); }
private double SecondDerivativeByFlux(SinglePhaseField field, int jCell, NodeSet nodeSet) { if (this.gradientX == null) { // Evaluate gradient gradientX = new SinglePhaseField(field.Basis, "gradientX"); gradientY = new SinglePhaseField(field.Basis, "gradientY"); gradientX.DerivativeByFlux(1.0, field, d: 0); gradientY.DerivativeByFlux(1.0, field, d: 1); if (this.patchRecoveryGradient) { gradientX = ShockFindingExtensions.PatchRecovery(gradientX); gradientY = ShockFindingExtensions.PatchRecovery(gradientY); } } if (this.hessianXX == null) { // Evaluate Hessian matrix hessianXX = new SinglePhaseField(field.Basis, "hessianXX"); hessianXY = new SinglePhaseField(field.Basis, "hessianXY"); hessianYX = new SinglePhaseField(field.Basis, "hessianYX"); hessianYY = new SinglePhaseField(field.Basis, "hessianYY"); hessianXX.DerivativeByFlux(1.0, gradientX, d: 0); hessianXY.DerivativeByFlux(1.0, gradientX, d: 1); hessianYX.DerivativeByFlux(1.0, gradientY, d: 0); hessianYY.DerivativeByFlux(1.0, gradientY, d: 1); if (this.patchRecoveryHessian) { hessianXX = ShockFindingExtensions.PatchRecovery(hessianXX); hessianXY = ShockFindingExtensions.PatchRecovery(hessianXY); hessianYX = ShockFindingExtensions.PatchRecovery(hessianYX); hessianYY = ShockFindingExtensions.PatchRecovery(hessianYY); } } if (this.resultGradX == null) { resultGradX = MultidimensionalArray.Create(1, 1); resultGradY = MultidimensionalArray.Create(1, 1); } if (this.resultHessXX == null) { resultHessXX = MultidimensionalArray.Create(1, 1); resultHessXY = MultidimensionalArray.Create(1, 1); resultHessYX = MultidimensionalArray.Create(1, 1); resultHessYY = MultidimensionalArray.Create(1, 1); } gradientX.Evaluate(jCell, 1, nodeSet, resultGradX); gradientY.Evaluate(jCell, 1, nodeSet, resultGradY); hessianXX.Evaluate(jCell, 1, nodeSet, resultHessXX); hessianXY.Evaluate(jCell, 1, nodeSet, resultHessXY); hessianYX.Evaluate(jCell, 1, nodeSet, resultHessYX); hessianYY.Evaluate(jCell, 1, nodeSet, resultHessYY); // Compute second derivative along curve double g_alpha_alpha = 2 * ((resultHessXX[0, 0] * resultGradX[0, 0] + resultHessXY[0, 0] * resultGradY[0, 0]) * resultGradX[0, 0] + (resultHessYX[0, 0] * resultGradX[0, 0] + resultHessYY[0, 0] * resultGradY[0, 0]) * resultGradY[0, 0]); if (g_alpha_alpha == 0.0 || g_alpha_alpha.IsNaN()) { throw new NotSupportedException("Second derivative is zero"); } return(g_alpha_alpha); }
/// <summary> /// Solves the Reinitialisation-problem locally, i.e. cell by cell. /// </summary> /// <param name="Accepted"> /// </param> /// <param name="Recalc"></param> /// <param name="phi_accepted"> /// Input: level-set field values for the accepted cells/vallues in these cells are unchanged on exit. /// Output: level-set-field values for cells in <see cref="Recalc"/>. /// </param> /// <param name="_sign"></param> /// <param name="gradPhi"></param> void LocalSolve(CellMask Accepted, CellMask Recalc, SinglePhaseField phi_accepted, VectorField <SinglePhaseField> gradPhi, double _sign, SinglePhaseField DiffusionCoeff) { // loop over 'Recalc' - cells... foreach (int jCell in Recalc.ItemEnum) { double mini, maxi; this.Stpw_geomSolver.Start(); this.lsGeom.LocalSolve(jCell, Accepted.GetBitMaskWithExternal(), phi_accepted, _sign, out mini, out maxi); this.Stpw_geomSolver.Stop(); this.Stpw_reinitSolver.Start(); this.lsEllipt.LocalSolve(jCell, Accepted.GetBitMaskWithExternal(), phi_accepted, gradPhi); this.Stpw_reinitSolver.Stop(); //bool q = LocalSolve_Iterative(jCell, Accepted.GetBitMask(), phi_accepted, gradPhi); //if(!q) // Console.WriteLine("Local solver not converged."); } }
/// <summary> /// Extension velocity on un-cut cells. /// </summary> /// <param name="Phi">Input; a signed-distance field.</param> /// <param name="Domain">Cells in which the extension velocity should be computed.</param> /// <param name="cut">Cells in which a value for the extension property is given, must be /// disjoint from and share a boundary with <paramref name="Domain"/>. /// </param> /// <param name="ExtProperty"> /// Input/Output: on cells in <paramref name="cut"/>, valid values for the extension property, /// i.e. boundary values for the extension problem. /// On exit, the result of the marching algorithm is stored in the <paramref name="Domain"/>-cells. /// Multiple extension properties can be constructed at once. /// </param> /// <param name="ExtPropertyMin"> /// Input/Output: Helper array to check the minimum/maximum principle of the extension velocity problem /// - 1st index: extension property index /// - 2nd index: cell index /// On entry, the minimum values for cells in <paramref name="cut"/> must contain valid entries. /// </param> /// <param name="ExtPropertyMax"> /// Analogous to <paramref name="ExtPropertyMin"/>. /// </param> /// <param name="GradPhi"> /// Helper variable to compute the gradient values of <paramref name="Phi"/>. /// </param> public void ConstructExtension(SinglePhaseField Phi, CellMask Domain, CellMask cut, ConventionalDGField[] ExtProperty, double[][] ExtPropertyMin, double[][] ExtPropertyMax, VectorField <SinglePhaseField> GradPhi, int TimestepNo, bool plotMarchingSteps = false) // { using (new FuncTrace()) { Tracer.InstrumentationSwitch = false; // lots of tracing on calls acting on singe cells causes massive overhead (up to 5x slower). Stpw_total.Start(); // check args and init // =================== //ExtVelSolver extVelSlv = null; ExtVelSolver_Geometric extVelSlv = null; if (ExtProperty != null) { if (ExtProperty.Length != ExtPropertyMin.Length) { throw new ArgumentException(); } if (ExtProperty.Length != ExtPropertyMax.Length) { throw new ArgumentException(); } //extVelSlv = new ExtVelSolver(ExtProperty[0].Basis); extVelSlv = new ExtVelSolver_Geometric(ExtProperty[0].Basis); } BitArray Accepted_Mutuable = cut.GetBitMask().CloneAs(); int J = this.GridDat.Cells.Count; int D = this.GridDat.SpatialDimension; int N = this.LevelSetBasis.Length; int[] DomainCellIndices = Domain.ItemEnum.ToArray(); double[] PhiAvg = new double[DomainCellIndices.Length]; { int L = DomainCellIndices.Length; for (int jSub = 0; jSub < L; jSub++) { int jCell = DomainCellIndices[jSub]; PhiAvg[jSub] = Phi.GetMeanValue(jCell); } Array.Sort(PhiAvg, DomainCellIndices); } if (this.GridDat.MpiSize > 1) { throw new NotSupportedException("Currently not MPI parallel."); } int[] PosDomain, NegDomain; { int median = 0; for (; median < DomainCellIndices.Length; median++) { if (PhiAvg[median] >= 0) { break; } } NegDomain = new int[median]; for (int i = 0; i < median; i++) { NegDomain[i] = DomainCellIndices[median - i - 1]; } PosDomain = new int[DomainCellIndices.Length - median]; Array.Copy(DomainCellIndices, median, PosDomain, 0, PosDomain.Length); Debug.Assert(PosDomain.Length + NegDomain.Length == DomainCellIndices.Length); } if (plotMarchingSteps) { this.plotter.setup(new DGField[] { Phi, ExtProperty[0], ExtProperty[1] }, TimestepNo); this.plotter.plotstep(Accepted_Mutuable); } // perform marching... // =================== // marching loop.. for (int iMinusPlus = -1; iMinusPlus <= 1; iMinusPlus += 2) { double _sign = iMinusPlus; int[] _Domain; switch (iMinusPlus) { case -1: _Domain = NegDomain; break; case +1: _Domain = PosDomain; break; default: throw new Exception(); } for (int iSub = 0; iSub < _Domain.Length; iSub++) { CellMask Accepted = new CellMask(this.GridDat, Accepted_Mutuable); int jCellAccpt = _Domain[iSub]; //this.Stpw_gradientEval.Start(); //gradModule.GradientUpdate(jCellAccpt, Acceped_Mutuable, Phi, GradPhi); //this.Stpw_gradientEval.Stop(); // solve for the extension properties // ---------------------------------- if (ExtProperty != null) { int[] Neighb, dummy33; GridDat.GetCellNeighbours(jCellAccpt, GetCellNeighbours_Mode.ViaEdges, out Neighb, out dummy33); // solve for each component seperately for (int iComp = 0; iComp < ExtProperty.Length; iComp++) { ExtPropertyMax[iComp][jCellAccpt] = -double.MaxValue; ExtPropertyMin[iComp][jCellAccpt] = double.MaxValue; foreach (int jNeig in Neighb) { if (Accepted_Mutuable[jNeig]) { ExtPropertyMax[iComp][jCellAccpt] = Math.Max(ExtPropertyMax[iComp][jCellAccpt], ExtPropertyMax[iComp][jNeig]); ExtPropertyMin[iComp][jCellAccpt] = Math.Min(ExtPropertyMin[iComp][jCellAccpt], ExtPropertyMin[iComp][jNeig]); } } this.Stpw_extVelSolver.Start(); //extVelSlv.ExtVelSolve_Far(Phi, GradPhi, ExtProperty[iComp], ref ExtPropertyMin[iComp][jCellAccpt], ref ExtPropertyMax[iComp][jCellAccpt], jCellAccpt, Accepted, _sign); extVelSlv.ExtVelSolve_Geometric(Phi, ExtProperty[iComp], Accepted_Mutuable, jCellAccpt, _sign); this.Stpw_extVelSolver.Start(); } } Accepted_Mutuable[jCellAccpt] = true; if (plotMarchingSteps) { this.plotter.plotstep(Accepted_Mutuable); } } } Tracer.InstrumentationSwitch = true; Stpw_total.Stop(); } }
/// <summary> /// Calculates the Torque around the center of mass /// </summary> /// <param name="U"></param> /// <param name="P"></param> /// <param name="momentFittingVariant"></param> /// <param name="muA"></param> /// <param name="particleRadius"></param> /// <returns></returns> static public void GetCellValues(VectorField <XDGField> U, XDGField P, double muA, double particleRadius, SinglePhaseField P_atIB, SinglePhaseField gradU_atIB, SinglePhaseField gradUT_atIB) { var LsTrk = U[0].Basis.Tracker; int D = LsTrk.GridDat.SpatialDimension; var UA = U.Select(u => u.GetSpeciesShadowField("A")).ToArray(); if (D > 2) { throw new NotImplementedException("Currently only 2D cases supported"); } int RequiredOrder = U[0].Basis.Degree * 3 + 2; //if (RequiredOrder > agg.HMForder) // throw new ArgumentException(); Console.WriteLine("Cell values calculated by: {0}, order = {1}", LsTrk.CutCellQuadratureType, RequiredOrder); ConventionalDGField pA = null; double circumference = new double(); pA = P.GetSpeciesShadowField("A"); for (int n = 0; n < 4; n++) { ScalarFunctionEx ErrFunc_CellVal = delegate(int j0, int Len, NodeSet Ns, MultidimensionalArray result) { int K = result.GetLength(1); // No nof Nodes MultidimensionalArray Grad_UARes = MultidimensionalArray.Create(Len, K, D, D);; MultidimensionalArray pARes = MultidimensionalArray.Create(Len, K); // Evaluate tangential velocity to level-set surface var Normals = LsTrk.DataHistories[0].Current.GetLevelSetNormals(Ns, j0, Len); for (int i = 0; i < D; i++) { UA[i].EvaluateGradient(j0, Len, Ns, Grad_UARes.ExtractSubArrayShallow(-1, -1, i, -1)); } pA.Evaluate(j0, Len, Ns, pARes); for (int j = 0; j < Len; j++) { for (int k = 0; k < K; k++) { double acc = 0.0; double acc2 = 0.0; switch (n) { case 0: // Pressure part acc += pARes[j, k] * Normals[j, k, 0]; acc *= -Normals[j, k, 1] * particleRadius; acc2 += pARes[j, k] * Normals[j, k, 1]; acc2 *= Normals[j, k, 0] * particleRadius; result[j, k] = acc + acc2; break; case 1: // GradU part acc -= (1 * muA) * Grad_UARes[j, k, 0, 0] * Normals[j, k, 0]; // Attention was 2 times acc -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 1]; acc *= -Normals[j, k, 1] * particleRadius; acc2 -= (1 * muA) * Grad_UARes[j, k, 1, 1] * Normals[j, k, 1]; acc2 -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 0]; acc2 *= Normals[j, k, 0] * particleRadius; result[j, k] = acc + acc2; break; case 2: // GradU_T part acc -= (1 * muA) * Grad_UARes[j, k, 0, 0] * Normals[j, k, 0]; // Attention was 2 times acc -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 1]; acc *= -Normals[j, k, 1] * particleRadius; acc2 -= (1 * muA) * Grad_UARes[j, k, 1, 1] * Normals[j, k, 1]; // Attention was 2 times acc2 -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 0]; acc2 *= Normals[j, k, 0] * particleRadius; result[j, k] = acc + acc2; break; case 3: // Standardization with radians result[j, k] = 1; break; default: throw new NotImplementedException(); } } } }; var SchemeHelper = LsTrk.GetXDGSpaceMetrics(new[] { LsTrk.GetSpeciesId("A") }, RequiredOrder, 1).XQuadSchemeHelper; // new XQuadSchemeHelper(LsTrk, momentFittingVariant, ); CellQuadratureScheme cqs = SchemeHelper.GetLevelSetquadScheme(0, LsTrk.Regions.GetCutCellMask()); CellQuadrature.GetQuadrature(new int[] { 1 }, LsTrk.GridDat, cqs.Compile(LsTrk.GridDat, RequiredOrder), delegate(int i0, int Length, QuadRule QR, MultidimensionalArray EvalResult) { ErrFunc_CellVal(i0, Length, QR.Nodes, EvalResult.ExtractSubArrayShallow(-1, -1, 0)); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { for (int i = 0; i < Length; i++) { switch (n) { case 0: P_atIB.SetMeanValue(i0, ResultsOfIntegration[i, 0]); break; case 1: gradU_atIB.SetMeanValue(i0, ResultsOfIntegration[i, 0]); break; case 2: gradUT_atIB.SetMeanValue(i0, ResultsOfIntegration[i, 0]); break; case 3: circumference += ResultsOfIntegration[i, 0]; P_atIB.SetMeanValue(i0, P_atIB.GetMeanValue(i0) / ResultsOfIntegration[i, 0]); gradU_atIB.SetMeanValue(i0, gradU_atIB.GetMeanValue(i0) / ResultsOfIntegration[i, 0]); gradUT_atIB.SetMeanValue(i0, gradUT_atIB.GetMeanValue(i0) / ResultsOfIntegration[i, 0]); break; default: throw new NotImplementedException(); } } } ).Execute(); } Console.WriteLine("Circle circumference: " + circumference); }
public bool LocalSolve(int jCell, BitArray AcceptedMask, SinglePhaseField Phi, double _sign, out double Min, out double Max) { Debug.Assert(_sign.Abs() == 1); // find all accepted neighbors // =========================== var NeighCells = this.GridDat.GetCellNeighboursViaEdges(jCell); int iKref = this.GridDat.Cells.GetRefElementIndex(jCell); int NN = NeighCells.Length; var NeighCellsK = new Tuple <int, int, int> [NN]; int NNK = 0; for (int nn = 0; nn < NN; nn++) { if (AcceptedMask[NeighCells[nn].Item1] == true) { NeighCellsK[NNK] = NeighCells[nn]; NNK++; } } // evaluate accepted neighbors // ============================ Min = double.MaxValue; Max = double.MinValue; var TrafoIdx = this.GridDat.Edges.Edge2CellTrafoIndex; int K = this.PhiEvalBuffer[0].GetLength(1); // Nodes per edge for (int nnk = 0; nnk < NNK; nnk++) // loop over accepted neighbours { int jNC = NeighCellsK[nnk].Item1; int iEdg = NeighCellsK[nnk].Item2; int InOrOut = NeighCellsK[nnk].Item3; int iTrafo = TrafoIdx[iEdg, InOrOut]; this.PhiEvalBuffer[nnk].Clear(); Phi.Evaluate(jNC, 1, this.EdgeNodes.GetVolumeNodeSet(this.GridDat, iTrafo), this.PhiEvalBuffer[nnk]); this.GridDat.TransformLocal2Global(this.EdgeNodes.GetVolumeNodeSet(this.GridDat, iTrafo), this.CellNodesGlobal[nnk], jNC); Max = Math.Max(Max, PhiEvalBuffer[nnk].Max()); Min = Math.Min(Min, PhiEvalBuffer[nnk].Min()); } var _QuadNodesGlobal = this.QuadNodesGlobal[iKref]; this.GridDat.TransformLocal2Global(this.DaRuleS[iKref].Nodes, _QuadNodesGlobal, jCell); if (_sign > 0) { Max += this.GridDat.Cells.h_max[jCell]; } else { Min -= this.GridDat.Cells.h_max[jCell]; } // perform projection of geometric reinit // ====================================== // temp storage double[] Y = new double[2]; double[] X = new double[2]; double[] X1 = new double[2]; double[] X2 = new double[2]; // basis values at cell quadrature nodes var BasisValues = this.LevelSetBasis_Geometric.CellEval(this.DaRuleS[iKref].Nodes, jCell, 1).ExtractSubArrayShallow(0, -1, -1); int NoOfQn = BasisValues.GetLength(0); // number of quadrature nodes // result at quadrature nodes var PhiAtQuadNodes = MultidimensionalArray.Create(NoOfQn); for (int iQn = NoOfQn - 1; iQn >= 0; iQn--) // loop over all quadrature nodes { Y[0] = _QuadNodesGlobal[iQn, 0]; Y[1] = _QuadNodesGlobal[iQn, 1]; double _dist_min1 = double.MaxValue, _phi_min1 = double.NaN; int _nnk_Min1 = int.MinValue, _k_min1 = int.MinValue; double _dist_min2 = double.MaxValue, _phi_min2 = double.NaN; int _nnk_Min2 = int.MinValue, _k_min2 = int.MinValue; // find closest point: brute force approach for (int nnk = 0; nnk < NNK; nnk++) // loop over all edges with known values { double dist_min1 = double.MaxValue, phi_min1 = double.NaN; int nnk_Min1 = int.MinValue, k_min1 = int.MinValue; double dist_min2 = double.MaxValue, phi_min2 = double.NaN; int nnk_Min2 = int.MinValue, k_min2 = int.MinValue; for (int k = 0; k < K; k++) // loop over all nodes on this edge { X[0] = this.CellNodesGlobal[nnk][k, 0]; X[1] = this.CellNodesGlobal[nnk][k, 1]; double phi = this.PhiEvalBuffer[nnk][0, k]; phi *= _sign; double dist = GenericBlas.L2Dist(X, Y) + phi; bool NoBlock = true; if (dist < dist_min1) { nnk_Min2 = nnk_Min1; k_min2 = k_min1; dist_min2 = dist_min1; phi_min2 = phi_min1; dist_min1 = dist; nnk_Min1 = nnk; k_min1 = k; phi_min1 = phi; NoBlock = false; } if (dist >= dist_min1 && dist < dist_min2 && NoBlock) { dist_min2 = dist; nnk_Min2 = nnk; k_min2 = k; phi_min2 = phi; } } if (dist_min1 < _dist_min1) { _dist_min1 = dist_min1; _k_min1 = k_min1; _nnk_Min1 = nnk_Min1; _phi_min1 = phi_min1; _dist_min2 = dist_min2; _k_min2 = k_min2; _nnk_Min2 = nnk_Min2; _phi_min2 = phi_min2; } } { Debug.Assert(_nnk_Min1 == _nnk_Min2); Debug.Assert(_k_min1 != _k_min2); double PhiMin1 = this.PhiEvalBuffer[_nnk_Min1][0, _k_min1]; double PhiMin2 = this.PhiEvalBuffer[_nnk_Min2][0, _k_min2]; X1[0] = this.CellNodesGlobal[_nnk_Min1][_k_min1, 0]; X1[1] = this.CellNodesGlobal[_nnk_Min1][_k_min1, 1]; X2[0] = this.CellNodesGlobal[_nnk_Min2][_k_min2, 0]; X2[1] = this.CellNodesGlobal[_nnk_Min2][_k_min2, 1]; } _dist_min1 *= _sign; PhiAtQuadNodes[iQn] = _dist_min1 * this.DaRuleS[iKref].Weights[iQn]; } // finalize projection & return // ============================ if (this.GridDat.Cells.IsCellAffineLinear(jCell)) { int N = this.LevelSetBasis_Geometric.GetLength(jCell); int N2 = Phi.Basis.GetLength(jCell); MultidimensionalArray Phi_1 = MultidimensionalArray.Create(N); double scale = this.GridDat.Cells.JacobiDet[jCell]; Phi_1.Multiply(scale, BasisValues, PhiAtQuadNodes, 0.0, "m", "km", "k"); for (int n = 0; n < N; n++) { Phi.Coordinates[jCell, n] = Phi_1[n]; } for (int n = N; n < N2; n++) { Phi.Coordinates[jCell, n] = 0; } } else { throw new NotImplementedException("not implemented for curved cells"); } return(true); }
public void Init(MultigridOperator op) { int D = op.GridData.SpatialDimension; CodName = (new string[] { "mom0", "mom1" }); Params = ArrayTools.Cat( VariableNames.Velocity0Vector(D)); DomName = ArrayTools.Cat(VariableNames.VelocityVector(D)); LocalOp = new SpatialOperator(DomName, Params, CodName, (A, B, C) => 4); for (int d = 0; d < D; d++) { LocalOp.EquationComponents["mom" + d].Add(new LocalDiffusiveFlux() { m_component = d, dt = m_dt, muA = m_muA }); } LocalOp.Commit(); //LocalMatrix = op.MassMatrix.CloneAs().ToMsrMatrix(); //LocalMatrix.Clear(); //var U0 = new VectorField<SinglePhaseField>(op.BaseGridProblemMapping Take(D).Select(F => (SinglePhaseField)F).ToArray()); UnsetteledCoordinateMapping test = new UnsetteledCoordinateMapping(op.BaseGridProblemMapping.BasisS.GetSubVector(0, D)); var U0 = ((BoSSS.Foundation.CoordinateMapping)op.Mapping.ProblemMapping).Fields.GetSubVector(0, 2); var empty = new SinglePhaseField[D]; LocalMatrix = LocalOp.ComputeMatrix(test, empty, test, time: m_dt); Uidx = op.Mapping.ProblemMapping.GetSubvectorIndices(true, D.ForLoop(i => i)); Pidx = op.Mapping.ProblemMapping.GetSubvectorIndices(true, D); int Upart = Uidx.Length; int Ppart = Pidx.Length; ConvDiff = new MsrMatrix(Upart, Upart, 1, 1); pGrad = new MsrMatrix(Upart, Ppart, 1, 1); divVel = new MsrMatrix(Ppart, Upart, 1, 1); var VelocityMass = new MsrMatrix(Upart, Upart, 1, 1); var leftChangeBasesVel = new MsrMatrix(Upart, Upart, 1, 1); var rightChangeBasesVel = new MsrMatrix(Upart, Upart, 1, 1); op.MassMatrix.AccSubMatrixTo(1.0, VelocityMass, Uidx, default(int[]), Uidx, default(int[]), default(int[]), default(int[])); op.LeftChangeOfBasis.AccSubMatrixTo(1.0, leftChangeBasesVel, Uidx, default(int[]), Uidx, default(int[]), default(int[]), default(int[])); op.RightChangeOfBasis.AccSubMatrixTo(1.0, rightChangeBasesVel, Uidx, default(int[]), Uidx, default(int[]), default(int[]), default(int[])); var temp = MsrMatrix.Multiply(leftChangeBasesVel, LocalMatrix); LocalMatrix = MsrMatrix.Multiply(temp, rightChangeBasesVel); var M = op.OperatorMatrix; M.AccSubMatrixTo(1.0, ConvDiff, Uidx, default(int[]), Uidx, default(int[]), default(int[]), default(int[])); M.AccSubMatrixTo(1.0, pGrad, Uidx, default(int[]), Pidx, default(int[]), default(int[]), default(int[])); M.AccSubMatrixTo(1.0, divVel, Pidx, default(int[]), Uidx, default(int[]), default(int[]), default(int[])); LocalMatrix.SaveToTextFileSparse("LocalConvDiffMatrix"); ConvDiff.SaveToTextFileSparse("ConvDiff"); op.MassMatrix.SaveToTextFileSparse("MassMatrix"); VelocityMass.SaveToTextFileSparse("VelocityMass"); }
/// <summary> /// % /// </summary> /// <param name="output"></param> /// <param name="input">input DG field; unchanged on </param> public void Perform(SinglePhaseField output, SinglePhaseField input) { if (!output.Basis.Equals(this.m_bOutput)) { throw new ArgumentException("output basis mismatch"); } if (!input.Basis.Equals(this.m_bInput)) { throw new ArgumentException("output basis mismatch"); } var GridDat = this.m_bOutput.GridDat; int N = m_bOutput.Length; int Nin = m_bInput.Length; int J = GridDat.iLogicalCells.NoOfLocalUpdatedCells; diagnosis = new double[N * J]; double[] RHS = new double[N]; double[] f2 = new double[N]; double[] g1 = new double[N]; MultidimensionalArray Trf = MultidimensionalArray.Create(N, N); input.MPIExchange(); for (int jCell = 0; jCell < J; jCell++) { Debug.Assert((this.AggregateBasisTrafo[jCell] != null) == (this.Stencils[jCell] != null)); //FullMatrix invMassM = this.InvMassMatrix[jCell]; MultidimensionalArray ExPolMtx = this.AggregateBasisTrafo[jCell]; if (ExPolMtx != null) { int[] Stencil_jCells = this.Stencils[jCell]; int K = Stencil_jCells.Length; Debug.Assert(Stencil_jCells[0] == jCell); var coordIn = input.Coordinates; //for(int n = Math.Min(Nin, N) - 1; n >= 0; n--) { // RHS[n] = coordIn[jCell, n]; //} RHS.ClearEntries(); g1.ClearEntries(); var oldCoords = input.Coordinates.GetRow(jCell); for (int k = 0; k < K; k++) { int jNeigh = Stencil_jCells[k]; for (int n = Math.Min(input.Basis.Length, N) - 1; n >= 0; n--) { f2[n] = input.Coordinates[jNeigh, n]; } if (Nlim != null && Nlim[jNeigh] > 0) { for (int n = Nlim[jNeigh]; n < RHS.Length; n++) { f2[n] = 0; } } for (int l = 0; l < N; l++) { double acc = 0; for (int i = 0; i < N; i++) { acc += ExPolMtx[k, i, l] * f2[i]; } RHS[l] += acc; } } if (Nlim != null && Nlim[jCell] > 0) { for (int n = Nlim[jCell]; n < RHS.Length; n++) { RHS[n] = 0; } } for (int l = 0; l < N; l++) { diagnosis[jCell * N + l] = RHS[l]; } //invMassM.gemv(1.0, RHS, 0.0, g1); Trf.Clear(); Trf.Acc(1.0, ExPolMtx.ExtractSubArrayShallow(0, -1, -1)); //Trf.Solve(FulCoords, AggCoords); Trf.GEMV(1.0, RHS, 0.0, g1); if (Nlim != null && Nlim[jCell] > 0) { for (int n = Nlim[jCell]; n < RHS.Length; n++) { g1[n] = 0; } } if (Nlim == null) { output.Coordinates.SetRow(jCell, g1); } else { if ((Nlim[jCell] <= 0 && this.notchangeunlim)) { output.Coordinates.SetRow(jCell, oldCoords); } else { output.Coordinates.SetRow(jCell, g1); } } } } output.MPIExchange(); }
//static void Main(string[] args) { // Application<CNSControl>._Main(args, false, null, () => new ViscousShockProfile()); //} protected override void SetInitial() { //base.SetInitial(); WorkingSet.ProjectInitialValues(SpeciesMap, base.Control.InitialValues_Evaluators); string viscosityLaw = Control.ViscosityLaw.ToString(); if (true) { DGField mpiRank = new SinglePhaseField(new Basis(GridData, 0), "rank"); m_IOFields.Add(mpiRank); for (int j = 0; j < GridData.Cells.NoOfLocalUpdatedCells; j++) { mpiRank.SetMeanValue(j, DatabaseDriver.MyRank); } } // exakte Lösung für 1D-Testfall int p = Control.MomentumDegree; CellQuadratureScheme scheme = new CellQuadratureScheme(true); int order = WorkingSet.Momentum[0].Basis.Degree * 2 + 2; int noOfNodesPerCell = scheme.Compile(GridData, order).First().Rule.NoOfNodes; int offset = Grid.CellPartitioning.i0 * noOfNodesPerCell; long K = Grid.NumberOfCells; string pathToData; string initial; if (Grid.Description.Contains("Unsteady")) { pathToData = @"P:\ViscousTerms\NS_exact\data\M4_" + viscosityLaw + "_unsteady\\"; initial = "1"; } else { pathToData = @"P:\ViscousTerms\NS_exact\data\M4_" + viscosityLaw + "_steady\\"; initial = "0"; } ScalarFunction func = new ScalarFunction( (MultidimensionalArray arrIn, MultidimensionalArray arrOut) => { string filename = "m" + initial + "_" + K.ToString() + "_" + Control.MomentumDegree.ToString() + ".txt"; double[] output; using (System.IO.StreamReader sr = new System.IO.StreamReader(pathToData + filename)) { output = sr.ReadToEnd(). Split(','). Select((str) => double.Parse(str, System.Globalization.CultureInfo.InvariantCulture.NumberFormat)). Skip(offset). Take(Grid.CellPartitioning.LocalLength * noOfNodesPerCell). ToArray(); }; output.CopyTo(arrOut.Storage, 0); } ); WorkingSet.Momentum[0].ProjectField(1.0, func, scheme); p = Control.EnergyDegree; scheme = new CellQuadratureScheme(true); order = WorkingSet.Energy.Basis.Degree * 2 + 2; noOfNodesPerCell = scheme.Compile(GridData, order).First().Rule.NoOfNodes; offset = Grid.CellPartitioning.i0 * noOfNodesPerCell; func = new ScalarFunction( (MultidimensionalArray arrIn, MultidimensionalArray arrOut) => { string filename = "rhoE" + initial + "_" + K.ToString() + "_" + p.ToString() + ".txt"; double[] output; using (System.IO.StreamReader sr = new System.IO.StreamReader(pathToData + filename)) { output = sr.ReadToEnd(). Split(','). Select((str) => double.Parse(str, System.Globalization.CultureInfo.InvariantCulture.NumberFormat)). Skip(offset). Take(Grid.CellPartitioning.LocalLength * noOfNodesPerCell). ToArray(); }; output.CopyTo(arrOut.Storage, 0); } ); WorkingSet.Energy.ProjectField(func); p = Control.DensityDegree; scheme = new CellQuadratureScheme(true); order = WorkingSet.Density.Basis.Degree * 2 + 2; noOfNodesPerCell = scheme.Compile(GridData, order).First().Rule.NoOfNodes; offset = Grid.CellPartitioning.i0 * noOfNodesPerCell; func = new ScalarFunction( (MultidimensionalArray arrIn, MultidimensionalArray arrOut) => { string filename = "rho" + initial + "_" + K.ToString() + "_" + p.ToString() + ".txt"; double[] output; using (System.IO.StreamReader sr = new System.IO.StreamReader(pathToData + filename)) { output = sr.ReadToEnd(). Split(','). Select((str) => double.Parse(str, System.Globalization.CultureInfo.InvariantCulture.NumberFormat)). Skip(offset). Take(Grid.CellPartitioning.LocalLength * noOfNodesPerCell). ToArray(); }; output.CopyTo(arrOut.Storage, 0); } ); WorkingSet.Density.ProjectField(func); WorkingSet.UpdateDerivedVariables(this, CellMask.GetFullMask(GridData)); //string guidString = "00000000-0000-0000-0000-000000000000"; //switch (K) { // case 32: // guidString = "6725b9fc-d072-44cd-ae72-196e8a692735"; // break; // case 64: // guidString = "cb714aec-4b4a-4dae-9e70-32861daa195f"; // break; // case 128: // guidString = "678dc736-bbf6-4a1e-86eb-4f80b2354748"; // break; //} //Guid tsGuid = new Guid(guidString); //var db = this.GetDatabase(); //string fieldName = "rho"; //ITimestepInfo tsi = db.Controller.DBDriver.LoadTimestepInfo(tsGuid, null, db); //previousRho = (SinglePhaseField)db.Controller.DBDriver.LoadFields(tsi, GridData, new[] { fieldName }).Single().CloneAs(); //diffRho = WorkingSet.Density.CloneAs(); //diffRho.Clear(); //fieldName = "rhoE"; //previousRhoE = (SinglePhaseField)db.Controller.DBDriver.LoadFields(tsi, GridData, new[] { fieldName }).Single().CloneAs(); //diffRhoE = WorkingSet.Energy.CloneAs(); //diffRhoE.Clear(); //fieldName = "m0"; //previousM0 = (SinglePhaseField)db.Controller.DBDriver.LoadFields(tsi, GridData, new[] { fieldName }).Single().CloneAs(); //diffM0 = WorkingSet.Momentum[0].CloneAs(); //diffM0.Clear(); //diffRho.Identification = "diffRho"; //diffM0.Identification = "diffM0"; //diffRhoE.Identification = "diffRhoE"; //previousRho.Identification = "Rho_analy"; //previousM0.Identification = "M0_analy"; //previousRhoE.Identification = "RhoE_analy"; //m_IOFields.Add(diffM0); //m_IOFields.Add(diffRho); //m_IOFields.Add(diffRhoE); //m_IOFields.Add(previousM0); //m_IOFields.Add(previousRho); //m_IOFields.Add(previousRhoE); }
/// <summary> /// /// </summary> /// <param name="output">output</param> /// <param name="cm">mask on which the filtering should be performed</param> /// <param name="input">input</param> /// <param name="mode"></param> public static void FilterStencilProjection(SinglePhaseField output, CellMask cm, SinglePhaseField input, PatchRecoveryMode mode) { var GridDat = input.GridDat; if (!object.ReferenceEquals(output.GridDat, GridDat)) { throw new ArgumentException(); } if (!object.ReferenceEquals(input.GridDat, GridDat)) { throw new ArgumentException(); } if (cm == null) { cm = CellMask.GetFullMask(GridDat); } switch (mode) { case PatchRecoveryMode.L2_unrestrictedDom: case PatchRecoveryMode.L2_restrictedDom: { var Mask = cm.GetBitMaskWithExternal(); bool RestrictToCellMask = mode == PatchRecoveryMode.L2_restrictedDom; foreach (int jCell in cm.ItemEnum) { int[] NeighCells, dummy; GridDat.GetCellNeighbours(jCell, GetCellNeighbours_Mode.ViaVertices, out NeighCells, out dummy); if (RestrictToCellMask == true) { NeighCells = NeighCells.Where(j => Mask[j]).ToArray(); } FilterStencilProjection(output, jCell, NeighCells, input); } return; } case PatchRecoveryMode.ChebychevInteroplation: { int P = input.Basis.Degree * 3; double[] ChebyNodes = ChebyshevNodes(P); NodeSet Nds = new NodeSet(GridDat.iGeomCells.RefElements[0], P, 1); Nds.SetColumn(0, ChebyNodes); Nds.LockForever(); Polynomial[] Polynomials = GetNodalPolynomials(ChebyNodes); MultidimensionalArray PolyVals = MultidimensionalArray.Create(P, P); for (int p = 0; p < P; p++) { Polynomials[p].Evaluate(PolyVals.ExtractSubArrayShallow(p, -1), Nds); for (int i = 0; i < P; i++) { Debug.Assert(Math.Abs(((i == p) ? 1.0 : 0.0) - PolyVals[p, i]) <= 1.0e-8); } } foreach (int jCell in cm.ItemEnum) { AffineTrafo T; var NodalValues = NodalPatchRecovery(jCell, ChebyNodes, input, out T); /* * MultidimensionalArray Nodes2D = MultidimensionalArray.Create(P, P, 2); * for (int i = 0; i < P; i++) { * for (int j = 0; j < P; j++) { * Nodes2D[i, j, 0] = ChebyNodes[i]; * Nodes2D[i, j, 1] = ChebyNodes[j]; * } * } * var _Nodes2D = Nodes2D.ResizeShallow(P*P, 2); * var xNodes = _Nodes2D.ExtractSubArrayShallow(new int[] { 0, 0 }, new int[] { P*P - 1, 0 }); * var yNodes = _Nodes2D.ExtractSubArrayShallow(new int[] { 0, 1 }, new int[] { P*P - 1, 1 }); * int K = P*P; * * MultidimensionalArray xPolyVals = MultidimensionalArray.Create(P, K); * MultidimensionalArray yPolyVals = MultidimensionalArray.Create(P, K); * * for (int p = 0; p < P; p++) { * Polynomials[p].Evaluate(yPolyVals.ExtractSubArrayShallow(p, -1), yNodes); * Polynomials[p].Evaluate(xPolyVals.ExtractSubArrayShallow(p, -1), xNodes); * * } * * * double ERR = 0; * for (int k = 0; k < K; k++) { * * double x = xNodes[k, 0]; * double y = yNodes[k, 0]; * * double acc = 0; * double BasisHits = 0.0; * for (int i = 0; i < P; i++) { * for (int j = 0; j < P; j++) { * * bool isNode = (Math.Abs(x - ChebyNodes[i]) < 1.0e-8) && (Math.Abs(y - ChebyNodes[j]) < 1.0e-8); * double BasisSoll = isNode ? 1.0 : 0.0; * if (isNode) * Console.WriteLine(); * * double Basis_ij = xPolyVals[i, k]*yPolyVals[j, k]; * * Debug.Assert((Basis_ij - BasisSoll).Abs() < 1.0e-8); * * BasisHits += Math.Abs(Basis_ij); * * acc += Basis_ij*NodalValues[i, j]; * } * } * Debug.Assert(Math.Abs(BasisHits - 1.0) < 1.0e-8); * * * * * double soll = yNodes[k,0]; * Debug.Assert((soll - acc).Pow2() < 1.0e-7); * ERR = (soll - acc).Pow2(); * } * Console.WriteLine(ERR); */ output.ProjectField(1.0, delegate(int j0, int Len, NodeSet Nodes, MultidimensionalArray result) { var K = Nodes.NoOfNodes; MultidimensionalArray NT = T.Transform(Nodes); var xNodes = new NodeSet(null, NT.ExtractSubArrayShallow(new int[] { 0, 0 }, new int[] { K - 1, 0 })); var yNodes = new NodeSet(null, NT.ExtractSubArrayShallow(new int[] { 0, 1 }, new int[] { K - 1, 1 })); MultidimensionalArray xPolyVals = MultidimensionalArray.Create(P, K); MultidimensionalArray yPolyVals = MultidimensionalArray.Create(P, K); for (int p = 0; p < P; p++) { Polynomials[p].Evaluate(yPolyVals.ExtractSubArrayShallow(p, -1), yNodes); Polynomials[p].Evaluate(xPolyVals.ExtractSubArrayShallow(p, -1), xNodes); } for (int k = 0; k < K; k++) { double acc = 0; for (int i = 0; i < P; i++) { for (int j = 0; j < P; j++) { acc += xPolyVals[i, k] * yPolyVals[j, k] * NodalValues[i, j]; } } result[0, k] = acc; } }, new CellQuadratureScheme(true, new CellMask(input.GridDat, Chunk.GetSingleElementChunk(jCell)))); } return; } } }
/// <summary> /// Update Forces and Torque acting from fluid onto the particle /// </summary> /// <param name="U"></param> /// <param name="P"></param> /// <param name="LsTrk"></param> /// <param name="muA"></param> public void UpdateForcesAndTorque(VectorField <SinglePhaseField> U, SinglePhaseField P, LevelSetTracker LsTrk, double muA, double dt, double fluidDensity, bool NotFullyCoupled) { if (skipForceIntegration) { skipForceIntegration = false; return; } HydrodynamicForces[0][0] = 0; HydrodynamicForces[0][1] = 0; HydrodynamicTorque[0] = 0; int RequiredOrder = U[0].Basis.Degree * 3 + 2; Console.WriteLine("Forces coeff: {0}, order = {1}", LsTrk.CutCellQuadratureType, RequiredOrder); double[] Forces = new double[SpatialDim]; SinglePhaseField[] UA = U.ToArray(); ConventionalDGField pA = null; pA = P; if (IncludeTranslation) { for (int d = 0; d < SpatialDim; d++) { void ErrFunc(int CurrentCellID, int Length, NodeSet Ns, MultidimensionalArray result) { int NumberOfNodes = result.GetLength(1); MultidimensionalArray Grad_UARes = MultidimensionalArray.Create(Length, NumberOfNodes, SpatialDim, SpatialDim); MultidimensionalArray pARes = MultidimensionalArray.Create(Length, NumberOfNodes); var Normals = LsTrk.DataHistories[0].Current.GetLevelSetNormals(Ns, CurrentCellID, Length); for (int i = 0; i < SpatialDim; i++) { UA[i].EvaluateGradient(CurrentCellID, Length, Ns, Grad_UARes.ExtractSubArrayShallow(-1, -1, i, -1), 0, 1); } pA.Evaluate(CurrentCellID, Length, Ns, pARes); for (int j = 0; j < Length; j++) { for (int k = 0; k < NumberOfNodes; k++) { result[j, k] = ForceIntegration.CalculateStressTensor(Grad_UARes, pARes, Normals, muA, k, j, this.SpatialDim, d); } } } var SchemeHelper = LsTrk.GetXDGSpaceMetrics(new[] { LsTrk.GetSpeciesId("A") }, RequiredOrder, 1).XQuadSchemeHelper; CellQuadratureScheme cqs = SchemeHelper.GetLevelSetquadScheme(0, CutCells_P(LsTrk)); CellQuadrature.GetQuadrature(new int[] { 1 }, LsTrk.GridDat, cqs.Compile(LsTrk.GridDat, RequiredOrder), delegate(int i0, int Length, QuadRule QR, MultidimensionalArray EvalResult) { ErrFunc(i0, Length, QR.Nodes, EvalResult.ExtractSubArrayShallow(-1, -1, 0)); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { Forces[d] = ParticleAuxillary.ForceTorqueSummationWithNeumaierArray(Forces[d], ResultsOfIntegration, Length); } ).Execute(); } } double Torque = 0; if (IncludeRotation) { void ErrFunc2(int j0, int Len, NodeSet Ns, MultidimensionalArray result) { int K = result.GetLength(1); // No nof Nodes MultidimensionalArray Grad_UARes = MultidimensionalArray.Create(Len, K, SpatialDim, SpatialDim);; MultidimensionalArray pARes = MultidimensionalArray.Create(Len, K); // Evaluate tangential velocity to level-set surface var Normals = LsTrk.DataHistories[0].Current.GetLevelSetNormals(Ns, j0, Len); for (int i = 0; i < SpatialDim; i++) { UA[i].EvaluateGradient(j0, Len, Ns, Grad_UARes.ExtractSubArrayShallow(-1, -1, i, -1), 0, 1); } pA.Evaluate(j0, Len, Ns, pARes); for (int j = 0; j < Len; j++) { MultidimensionalArray tempArray = Ns.CloneAs(); LsTrk.GridDat.TransformLocal2Global(Ns, tempArray, j0 + j); for (int k = 0; k < K; k++) { result[j, k] = ForceIntegration.CalculateTorqueFromStressTensor2D(Grad_UARes, pARes, Normals, tempArray, muA, k, j, Position[0]); } } } var SchemeHelper2 = LsTrk.GetXDGSpaceMetrics(new[] { LsTrk.GetSpeciesId("A") }, RequiredOrder, 1).XQuadSchemeHelper; CellQuadratureScheme cqs2 = SchemeHelper2.GetLevelSetquadScheme(0, CutCells_P(LsTrk)); CellQuadrature.GetQuadrature(new int[] { 1 }, LsTrk.GridDat, cqs2.Compile(LsTrk.GridDat, RequiredOrder), delegate(int i0, int Length, QuadRule QR, MultidimensionalArray EvalResult) { ErrFunc2(i0, Length, QR.Nodes, EvalResult.ExtractSubArrayShallow(-1, -1, 0)); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { Torque = ParticleAuxillary.ForceTorqueSummationWithNeumaierArray(Torque, ResultsOfIntegration, Length); } ).Execute(); } // add gravity { Forces[1] += (particleDensity - fluidDensity) * Area_P * GravityVertical; } // Sum forces and moments over all MPI processors // ============================================== { int NoOfVars = 1 + SpatialDim; double[] StateBuffer = new double[NoOfVars]; StateBuffer[0] = Torque; for (int d = 0; d < SpatialDim; d++) { StateBuffer[1 + d] = Forces[d]; } double[] GlobalStateBuffer = StateBuffer.MPISum(); Torque = GlobalStateBuffer[0]; for (int d = 0; d < SpatialDim; d++) { Forces[d] = GlobalStateBuffer[1 + d]; } } if (neglectAddedDamping == false) { double fest = Forces[0]; Forces[0] = Forces[0] + AddedDampingCoefficient * dt * (AddedDampingTensor[0, 0] * TranslationalAcceleration[0][0] + AddedDampingTensor[1, 0] * TranslationalAcceleration[0][1] + AddedDampingTensor[0, 2] * RotationalAcceleration[0]); Forces[1] = Forces[1] + AddedDampingCoefficient * dt * (AddedDampingTensor[0, 1] * TranslationalAcceleration[0][0] + AddedDampingTensor[1, 1] * TranslationalAcceleration[0][1] + AddedDampingTensor[1, 2] * RotationalAcceleration[0]); Torque += AddedDampingCoefficient * dt * (AddedDampingTensor[2, 0] * TranslationalAcceleration[0][0] + AddedDampingTensor[2, 1] * TranslationalAcceleration[0][1] + AddedDampingTensor[2, 2] * RotationalAcceleration[0]); } if (iteration_counter_P == -1 || NotFullyCoupled || iteration_counter_P == 250 || stupidcounter == 0) { Console.WriteLine(); if (iteration_counter_P == 1) { Console.WriteLine("First iteration of the current timestep, all relaxation factors are set to 1"); } if (iteration_counter_P == 250) { Console.WriteLine("250 iterations, I'm trying to jump closer to the real solution"); } for (int d = 0; d < SpatialDim; d++) { HydrodynamicForces[0][d] = 0; if (Math.Abs(Forces[d]) < ForceAndTorque_convergence * 1e-2 && ClearSmallValues == true) { Forces[d] = 0; } HydrodynamicForces[0][d] = Forces[d]; } HydrodynamicTorque[0] = 0; if (Math.Abs(Torque) < ForceAndTorque_convergence * 1e-2 && ClearSmallValues == true) { Torque = 0; } HydrodynamicTorque[0] = Torque; stupidcounter = 1; } else { double[] RelaxatedForceAndTorque = Underrelaxation.RelaxatedForcesAndTorque(Forces, Torque, ForcesPrevIteration, TorquePrevIteration, ForceAndTorque_convergence, underrelaxation_factor, ClearSmallValues, AddaptiveUnderrelaxation, AverageDistance, iteration_counter_P); for (int d = 0; d < SpatialDim; d++) { HydrodynamicForces[0][d] = RelaxatedForceAndTorque[d]; } HydrodynamicTorque[0] = RelaxatedForceAndTorque[SpatialDim]; } //for (int d = 0; d < SpatialDim; d++)// changes sign depending on the sign of Forces[d], should increase the convergence rate. (testing needed) //{ // if (Math.Abs(HydrodynamicForces[0][d] - Forces[0]) > Math.Abs(Forces[d])) // { // HydrodynamicForces[0][d] *= -1; // } //} if (double.IsNaN(HydrodynamicForces[0][0]) || double.IsInfinity(HydrodynamicForces[0][0])) { throw new ArithmeticException("Error trying to calculate hydrodynamic forces (x). Value: " + HydrodynamicForces[0][0]); } if (double.IsNaN(HydrodynamicForces[0][1]) || double.IsInfinity(HydrodynamicForces[0][1])) { throw new ArithmeticException("Error trying to calculate hydrodynamic forces (y). Value: " + HydrodynamicForces[0][1]); } if (double.IsNaN(HydrodynamicTorque[0]) || double.IsInfinity(HydrodynamicTorque[0])) { throw new ArithmeticException("Error trying to calculate hydrodynamic torque. Value: " + HydrodynamicTorque[0]); } }
/// <summary> /// /// </summary> /// <param name="jCell"></param> /// <param name="Stencil_jCell"></param> /// <param name="input">input</param> /// <param name="output">result</param> static void FilterStencilProjection(SinglePhaseField output, int jCell, int[] Stencil_jCells, ConventionalDGField input) { // implementation notes: Fk, persönliche Notizen, 17oct13 // ------------------------------------------------------ if (output.Basis.Degree < output.Basis.Degree) { throw new ArgumentException(); } int N = output.Basis.Length; // Basis dimension of 'g' in cell 'jCell' int K = Stencil_jCells.Length; // number of cells in stencil var ExPolMtx = MultidimensionalArray.Create(K, N, N); int[,] CellPairs = new int[K, 2]; for (int k = 0; k < K; k++) { CellPairs[k, 0] = jCell; CellPairs[k, 1] = Stencil_jCells[k]; } output.Basis.GetExtrapolationMatrices(CellPairs, ExPolMtx, null); MultidimensionalArray MassMatrix = MultidimensionalArray.Create(N, N); MassMatrix.AccEye(1.0); // Mass matrix in jCell itself //for (int k = 0; k < K; k++) { // over stencil members ... // for (int l = 0; l < N; l++) { // over rows of mass matrix ... // for (int m = 0; m < N; m++) { // over columns of mass matrix ... // double mass_lm = 0.0; // for (int i = 0; i < N; i++) { // mass_lm += ExPolMtx[k, i, m]*ExPolMtx[k, i, l]; // } // MassMatrix[l, m] += mass_lm; // } // } //} MassMatrix.Multiply(1.0, ExPolMtx, ExPolMtx, 1.0, "lm", "kim", "kil"); double[] RHS = new double[N]; double[] f2 = new double[N]; for (int n = Math.Min(input.Basis.Length, N) - 1; n >= 0; n--) { RHS[n] = input.Coordinates[jCell, n]; } for (int k = 0; k < K; k++) { for (int n = Math.Min(input.Basis.Length, N) - 1; n >= 0; n--) { f2[n] = input.Coordinates[Stencil_jCells[k], n]; } for (int l = 0; l < N; l++) { double acc = 0; for (int i = 0; i < N; i++) { acc += ExPolMtx[k, i, l] * f2[i]; } RHS[l] += acc; } } double[] g1 = new double[N]; MassMatrix.Solve(g1, RHS); output.Coordinates.SetRow(jCell, g1); }
/// <summary> /// Deprecated utility, writes matrices for spectral element method. /// </summary> public void WriteSEMMatrices() { int kSEM; // nodes per edge if (this.Grid.RefElements[0] is Triangle) { kSEM = this.T.Basis.Degree + 1; } else if (this.Grid.RefElements[0] is Square) { switch (this.T.Basis.Degree) { case 2: kSEM = 2; break; case 3: kSEM = 2; break; case 4: kSEM = 3; break; case 5: kSEM = 3; break; case 6: kSEM = 4; break; default: throw new NotSupportedException(); } } else { throw new NotImplementedException(); } SpecFemBasis SEM_basis = new SpecFemBasis((GridData)(this.GridData), kSEM); SEM_basis.CellNodes[0].SaveToTextFile("NODES_SEM" + kSEM + ".txt"); SEM_basis.MassMatrix.SaveToTextFileSparse("MASS_SEM" + kSEM + ".txt"); var Mod2Nod = SEM_basis.GetModal2NodalOperator(this.T.Basis); var Nod2Mod = SEM_basis.GetNodal2ModalOperator(this.T.Basis); Mod2Nod.SaveToTextFileSparse("MODAL" + this.T.Basis.Degree + "_TO_NODAL" + kSEM + ".txt"); Nod2Mod.SaveToTextFileSparse("NODAL" + kSEM + "_TO_MODAL" + this.T.Basis.Degree + ".txt"); this.LaplaceMtx.SaveToTextFileSparse("OPERATOR" + this.T.Basis.Degree + ".txt"); { var TEST = this.T.CloneAs(); TEST.Clear(); TEST.ProjectField((_2D)((x, y) => x * y)); int J = this.GridData.iGeomCells.Count; int N = SEM_basis.NodesPerCell[0]; MultidimensionalArray TEST_at_NODES = MultidimensionalArray.Create(J, N); TEST.Evaluate(0, J, SEM_basis.CellNodes[0], TEST_at_NODES); double[] CompareAtNodes = new double[J * N]; Mod2Nod.SpMVpara(1.0, TEST.CoordinateVector, 0.0, CompareAtNodes); double[] gretchen = TEST_at_NODES.ResizeShallow(J * N).To1DArray(); double fdist = GenericBlas.L2Dist(gretchen, CompareAtNodes); Debug.Assert(fdist < 1.0e-8); double[] bak = new double[TEST.Mapping.LocalLength]; Nod2Mod.SpMVpara(1.0, CompareAtNodes, 0.0, bak); double hdist = GenericBlas.L2Dist(bak, TEST.CoordinateVector); Debug.Assert(hdist < 1.0e-8); } var Scatter = SEM_basis.GetNodeScatterMatrix(); Scatter.SaveToTextFileSparse("SCATTER_SEM" + kSEM + ".txt"); { var SEMTEST = new SpecFemField(SEM_basis); var csem = SEMTEST.Coordinates; Random r = new Random(666); for (int i = 0; i < csem.GetLength(0); i++) { csem[i] = r.NextDouble(); } var DG_TEST0 = new SinglePhaseField(SEMTEST.Basis.ContainingDGBasis); var DG_TEST = this.T.CloneAs(); DG_TEST.Clear(); SEMTEST.AccToDGField(1.0, DG_TEST0); DG_TEST.AccLaidBack(1.0, DG_TEST0); double[] S2 = new double[Scatter.RowPartitioning.LocalLength]; double[] S3 = new double[DG_TEST.Mapping.LocalLength]; Scatter.SpMVpara(1.0, csem.To1DArray(), 0.0, S2); Nod2Mod.SpMVpara(1.0, S2, 0.0, S3); double gdist = GenericBlas.L2Dist(S3, DG_TEST.CoordinateVector); Debug.Assert(gdist < 1.0e-8); } }
/// <summary> /// /// </summary> /// <param name="jCell"></param> /// <param name="AcceptedMask"></param> /// <param name="Phi"></param> /// <param name="gradPhi"></param> /// <param name="__DiffusionCoeff">Output: if artificial diffusion is turned</param> /// <param name="MaxAllowedPhi">Input: upper threshold for the values of <paramref name="Phi"/> in cell <see cref="jCell"/>.</param> /// <param name="MinAllowedPhi">Input: lower threshold for the values of <paramref name="Phi"/> in cell <see cref="jCell"/>.</param> /// <returns></returns> public bool LocalSolve_Iterative(int jCell, BitArray AcceptedMask, SinglePhaseField Phi, VectorField <SinglePhaseField> gradPhi, SinglePhaseField __DiffusionCoeff, double MaxAllowedPhi, double MinAllowedPhi) { //this.LocalSolve_Geometric(jCell, AcceptedMask, Phi, +1, out MinAllowedPhi, out MaxAllowedPhi);) { int N = this.LevelSetBasis.GetLength(jCell); int i0G = this.LevelSetMapping.GlobalUniqueCoordinateIndex(0, jCell, 0); int i0L = this.LevelSetMapping.LocalUniqueCoordinateIndex(0, jCell, 0); SinglePhaseField __AcceptedMask = new SinglePhaseField(new Basis(this.GridDat, 0), "accepted"); for (int j = 0; j < AcceptedMask.Length; j++) { __AcceptedMask.SetMeanValue(j, AcceptedMask[j] ? 1.0 : 0.0); } // subgrid on which we are working, consisting only of one cell SubGrid jCellGrid = new SubGrid(new CellMask(this.GridDat, Chunk.GetSingleElementChunk(jCell))); // create spatial operator SpatialOperator.Evaluator evo; { SpatialOperator op = new SpatialOperator(1, 2, 1, QuadOrderFunc.NonLinear(2), "Phi", "dPhi_dx0", "dPhi_dx1", "cod1"); op.EquationComponents["cod1"].Add(new ReinitOperator()); op.Commit(); evo = op.GetEvaluatorEx(Phi.Mapping.Fields, gradPhi.Mapping.Fields, Phi.Mapping, edgeQrCtx: (new EdgeQuadratureScheme(domain: EdgeMask.GetEmptyMask(this.GridDat))), volQrCtx: (new CellQuadratureScheme(domain: jCellGrid.VolumeMask)), subGridBoundaryTreatment: SpatialOperator.SubGridBoundaryModes.InnerEdge); } // create artificial diffusion operator MultidimensionalArray DiffMtx; double[] DiffRhs; SpatialOperator dop; { double penaltyBase = this.LevelSetBasis.Degree + 2; penaltyBase = penaltyBase.Pow2(); dop = (new ArtificialViscosity(AcceptedMask, penaltyBase, GridDat.Cells.h_min, jCell, -1.0)).Operator(1); MsrMatrix _DiffMtx = new MsrMatrix(this.LevelSetMapping, this.LevelSetMapping); double[] _DiffRhs = new double[this.LevelSetMapping.LocalLength]; dop.ComputeMatrixEx(this.LevelSetMapping, new DGField[] { Phi, null, null }, this.LevelSetMapping, _DiffMtx, _DiffRhs, OnlyAffine: false, edgeQuadScheme: (new EdgeQuadratureScheme(domain: jCellGrid.AllEdgesMask)), volQuadScheme: (new CellQuadratureScheme(domain: jCellGrid.VolumeMask))); // extract matrix for 'jCell' DiffMtx = MultidimensionalArray.Create(N, N); DiffRhs = new double[N]; for (int n = 0; n < N; n++) { #if DEBUG int Lr; int[] row_cols = null; double[] row_vals = null; Lr = _DiffMtx.GetRow(i0G + n, ref row_cols, ref row_vals); for (int lr = 0; lr < Lr; lr++) { int ColIndex = row_cols[lr]; double Value = row_vals[lr]; Debug.Assert((ColIndex >= i0G && ColIndex < i0G + N) || (Value == 0.0), "Matrix is expected to be block-diagonal."); } #endif for (int m = 0; m < N; m++) { DiffMtx[n, m] = _DiffMtx[i0G + n, i0G + m]; } DiffRhs[n] = _DiffRhs[i0L + n]; } #if DEBUG var Test = DiffMtx.Transpose(); Test.Acc(-1.0, DiffMtx); Debug.Assert(Test.InfNorm() <= 1.0e-8); #endif } // find 'good' initial value by geometric solve AND // thresholds for the maximum an minimal value of Phi double Range = MaxAllowedPhi - MinAllowedPhi; MinAllowedPhi -= 0.1 * Range; MaxAllowedPhi += 0.1 * Range; // timestep for pseudo-timestepping double dt = 0.5 * this.GridDat.Cells.h_min[jCell] / (((double)(this.LevelSetBasis.Degree)).Pow2()); DGField[] PlotFields = new DGField[] { Phi, gradPhi[0], gradPhi[1], __DiffusionCoeff, __AcceptedMask }; //Tecplot.Tecplot.PlotFields(Params, "itt_0", "EllipicReinit", 0, 3); double[] PrevVal = new double[N]; double[] NextVal = new double[N]; Phi.Coordinates.GetRow(jCell, PrevVal); // pseudo-timestepping //if(jCell == 80) // Tecplot.Tecplot.PlotFields(PlotFields, "itt_0", "EllipicReinit", 0, 3); //Console.Write(" Local solve cell " + jCell + " ... "); bool converged = false; double DiffusionCoeff = 0; int IterGrowCount = 0; // number of iterations in which the residual grew double LastResi = double.NaN; for (int iIter = 0; iIter < 1000; iIter++) { //Console.Write(" Local solve iteration " + iIter + " ... "); PerformRKstep(dt, jCell, AcceptedMask, Phi, gradPhi, evo); __DiffusionCoeff.SetMeanValue(jCell, DiffusionCoeff); if (jCell == 80) { DiffusionCoeff = 0.1; } if (DiffusionCoeff > 0) { //Console.WriteLine(" Diffusion on."); double[] _DiffRhs = new double[this.LevelSetMapping.LocalLength]; dop.ComputeMatrixEx(this.LevelSetMapping, new DGField[] { Phi, gradPhi[0], gradPhi[1] }, this.LevelSetMapping, default(MsrMatrix), _DiffRhs, OnlyAffine: true, edgeQuadScheme: (new EdgeQuadratureScheme(domain: jCellGrid.AllEdgesMask)), volQuadScheme: (new CellQuadratureScheme(domain: CellMask.GetEmptyMask(this.GridDat)))); // extract matrix for 'jCell' for (int n = 0; n < N; n++) { DiffRhs[n] = _DiffRhs[i0L + n]; } PerformArtificialDiffusion(dt * DiffusionCoeff, jCell, Phi, DiffMtx, DiffRhs); } Phi.Coordinates.GetRow(jCell, NextVal); double resi = Math.Sqrt(GenericBlas.L2DistPow2(NextVal, PrevVal) / GenericBlas.L2NormPow2(PrevVal)); double[] A = NextVal; NextVal = PrevVal; PrevVal = A; if (iIter > 0 && resi > LastResi) { IterGrowCount++; } else { IterGrowCount = 0; } LastResi = resi; if (resi < 1.0e-10) { converged = true; break; } double maxPhi, minPhi; Phi.GetExtremalValuesInCell(out minPhi, out maxPhi, jCell); bool MinAlarm = minPhi < MinAllowedPhi; bool Maxalarm = maxPhi > MaxAllowedPhi; bool GrowAlarm = IterGrowCount > 4; bool IterAlarm = iIter >= 50; if (MinAlarm || Maxalarm || GrowAlarm) { // Diffusion coefficient should be increased if (DiffusionCoeff == 0) { DiffusionCoeff = 1.0e-2; } else { if (DiffusionCoeff < 1.0e3) { DiffusionCoeff *= 2; } } //Console.WriteLine(" increasing Diffusion: {0}, Alarms : {1}{2}{3}{4}", DiffusionCoeff, MinAlarm ? 1 : 0, Maxalarm ? 1 : 0, GrowAlarm ? 1 : 0, IterAlarm ? 1 : 0); } //if(jCell == 80 && iIter < 100) // Tecplot.Tecplot.PlotFields(PlotFields, "itt_" + (iIter + 1), "EllipicReinit", iIter + 1, 3); } return(converged); }
/// <summary> /// Obtaining the time integrated spatial discretization of the reinitialization equation in a narrow band around the zero level set, based on a Godunov's numerical Hamiltonian calculation /// </summary> /// <param name="LS"> The level set function </param> /// <param name="Restriction"> The narrow band around the zero level set </param> /// <param name="NumberOfTimesteps"> /// maximum number of pseudo-timesteps /// </param> /// <param name="thickness"> /// The smoothing width of the signum function. /// This is the main stabilization parameter for re-initialization. /// It should be set to approximately 3 cells. /// </param> /// <param name="TimestepSize"> /// size of the pseudo-timestep /// </param> public void ReInitialize(LevelSet LS, SubGrid Restriction, double thickness, double TimestepSize, int NumberOfTimesteps) { using (var tr = new FuncTrace()) { // log parameters: tr.Info("thickness: " + thickness.ToString(NumberFormatInfo.InvariantInfo)); tr.Info("TimestepSize: " + TimestepSize.ToString(NumberFormatInfo.InvariantInfo)); tr.Info("NumberOfTimesteps: " + NumberOfTimesteps); ExplicitEuler TimeIntegrator; SpatialOperator SO; Func <int[], int[], int[], int> QuadratureOrder = QuadOrderFunc.NonLinear(3); if (m_ctx.SpatialDimension == 2) { SO = new SpatialOperator(1, 5, 1, QuadratureOrder, new string[] { "LS", "LSCGV", "LSDG[0]", "LSUG[0]", "LSDG[1]", "LSUG[1]", "Result" }); SO.EquationComponents["Result"].Add(new GodunovHamiltonian(m_ctx, thickness)); SO.Commit(); TimeIntegrator = new RungeKutta(m_Scheme, SO, new CoordinateMapping(LS), new CoordinateMapping(LSCGV, LSDG[0], LSUG[0], LSDG[1], LSUG[1]), sgrd: Restriction); } else { SO = new SpatialOperator(1, 7, 1, QuadratureOrder, new string[] { "LS", "LSCGV", "LSDG[0]", "LSUG[0]", "LSDG[1]", "LSUG[1]", "LSDG[2]", "LSUG[2]", "Result" }); SO.EquationComponents["Result"].Add(new GodunovHamiltonian(m_ctx, thickness)); SO.Commit(); TimeIntegrator = new RungeKutta(m_Scheme, SO, new CoordinateMapping(LS), new CoordinateMapping(LSCGV, LSDG[0], LSUG[0], LSDG[1], LSUG[1], LSDG[2], LSUG[2]), sgrd: Restriction); } // Calculating the gradients in each sub-stage of a Runge-Kutta integration procedure ExplicitEuler.ChangeRateCallback EvalGradients = delegate(double t1, double t2) { LSUG.Clear(); CalculateLevelSetGradient(LS, LSUG, "Upwind", Restriction); LSDG.Clear(); CalculateLevelSetGradient(LS, LSDG, "Downwind", Restriction); LSCG.Clear(); CalculateLevelSetGradient(LS, LSCG, "Central", Restriction); LSCGV.Clear(); var VolMask = (Restriction != null) ? Restriction.VolumeMask : null; LSCGV.ProjectAbs(1.0, VolMask, LSCG.ToArray()); }; TimeIntegrator.OnBeforeComputeChangeRate += EvalGradients; { EvalGradients(0, 0); var GodunovResi = new SinglePhaseField(LS.Basis, "Residual"); SO.Evaluate(1.0, 0.0, LS.Mapping, TimeIntegrator.ParameterMapping.Fields, GodunovResi.Mapping, Restriction); //Tecplot.Tecplot.PlotFields(ArrayTools.Cat<DGField>( LSUG, LSDG, LS, GodunovResi), "Residual", 0, 3); } // pseudo-timestepping // =================== double factor = 1.0; double time = 0; LevelSet prevLevSet = new LevelSet(LS.Basis, "prevLevSet"); CellMask RestrictionMask = (Restriction == null) ? null : Restriction.VolumeMask; for (int i = 0; (i < NumberOfTimesteps); i++) { tr.Info("Level set reinitialization pseudo-timestepping, timestep " + i); // backup old Levelset // ------------------- prevLevSet.Clear(); prevLevSet.Acc(1.0, LS, RestrictionMask); // time integration // ---------------- double dt = TimestepSize * factor; tr.Info("dt = " + dt.ToString(NumberFormatInfo.InvariantInfo) + " (factor = " + factor.ToString(NumberFormatInfo.InvariantInfo) + ")"); TimeIntegrator.Perform(dt); time += dt; // change norm // ------ prevLevSet.Acc(-1.0, LS, RestrictionMask); double ChangeNorm = prevLevSet.L2Norm(RestrictionMask); Console.WriteLine("Reinit: PseudoTime: {0} - Changenorm: {1}", i, ChangeNorm); //Tecplot.Tecplot.PlotFields(new SinglePhaseField[] { LS }, m_ctx, "Reinit-" + i, "Reinit-" + i, i, 3); } //*/ } }
/// <summary> /// Calculates the drag (x-component) and lift (y-component) forces acting on a IB contour /// </summary> /// <param name="U"></param> /// <param name="P"></param> /// <param name="muA"></param> /// <returns></returns> static public double[] GetForces(VectorField <SinglePhaseField> U, SinglePhaseField P, LevelSetTracker LsTrk, double muA) { int D = LsTrk.GridDat.SpatialDimension; // var UA = U.Select(u => u.GetSpeciesShadowField("A")).ToArray(); var UA = U.ToArray(); int RequiredOrder = U[0].Basis.Degree * 3 + 2; //int RequiredOrder = LsTrk.GetXQuadFactoryHelper(momentFittingVariant).GetCachedSurfaceOrders(0).Max(); //Console.WriteLine("Order reduction: {0} -> {1}", _RequiredOrder, RequiredOrder); //if (RequiredOrder > agg.HMForder) // throw new ArgumentException(); Console.WriteLine("Forces coeff: {0}, order = {1}", LsTrk.CutCellQuadratureType, RequiredOrder); ConventionalDGField pA = null; //pA = P.GetSpeciesShadowField("A"); pA = P; double[] forces = new double[D]; for (int d = 0; d < D; d++) { ScalarFunctionEx ErrFunc = delegate(int j0, int Len, NodeSet Ns, MultidimensionalArray result) { int K = result.GetLength(1); // No nof Nodes MultidimensionalArray Grad_UARes = MultidimensionalArray.Create(Len, K, D, D);; MultidimensionalArray pARes = MultidimensionalArray.Create(Len, K); // Evaluate tangential velocity to level-set surface var Normals = LsTrk.DataHistories[0].Current.GetLevelSetNormals(Ns, j0, Len); for (int i = 0; i < D; i++) { UA[i].EvaluateGradient(j0, Len, Ns, Grad_UARes.ExtractSubArrayShallow(-1, -1, i, -1), 0, 1); } pA.Evaluate(j0, Len, Ns, pARes); if (LsTrk.GridDat.SpatialDimension == 2) { for (int j = 0; j < Len; j++) { for (int k = 0; k < K; k++) { double acc = 0.0; // pressure switch (d) { case 0: acc += pARes[j, k] * Normals[j, k, 0]; acc -= (2 * muA) * Grad_UARes[j, k, 0, 0] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 1]; break; case 1: acc += pARes[j, k] * Normals[j, k, 1]; acc -= (2 * muA) * Grad_UARes[j, k, 1, 1] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 0]; break; default: throw new NotImplementedException(); } result[j, k] = acc; } } } else { for (int j = 0; j < Len; j++) { for (int k = 0; k < K; k++) { double acc = 0.0; // pressure switch (d) { case 0: acc += pARes[j, k] * Normals[j, k, 0]; acc -= (2 * muA) * Grad_UARes[j, k, 0, 0] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 0, 2] * Normals[j, k, 2]; acc -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 2, 0] * Normals[j, k, 2]; break; case 1: acc += pARes[j, k] * Normals[j, k, 1]; acc -= (2 * muA) * Grad_UARes[j, k, 1, 1] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 1, 2] * Normals[j, k, 2]; acc -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 2, 1] * Normals[j, k, 2]; break; case 2: acc += pARes[j, k] * Normals[j, k, 2]; acc -= (2 * muA) * Grad_UARes[j, k, 2, 2] * Normals[j, k, 2]; acc -= (muA) * Grad_UARes[j, k, 2, 0] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 2, 1] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 0, 2] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 1, 2] * Normals[j, k, 1]; break; default: throw new NotImplementedException(); } result[j, k] = acc; } } } }; var SchemeHelper = LsTrk.GetXDGSpaceMetrics(new[] { LsTrk.GetSpeciesId("A") }, RequiredOrder, 1).XQuadSchemeHelper; //var SchemeHelper = new XQuadSchemeHelper(LsTrk, momentFittingVariant, ); CellQuadratureScheme cqs = SchemeHelper.GetLevelSetquadScheme(0, LsTrk.Regions.GetCutCellMask()); CellQuadrature.GetQuadrature(new int[] { 1 }, LsTrk.GridDat, cqs.Compile(LsTrk.GridDat, RequiredOrder), // agg.HMForder), delegate(int i0, int Length, QuadRule QR, MultidimensionalArray EvalResult) { ErrFunc(i0, Length, QR.Nodes, EvalResult.ExtractSubArrayShallow(-1, -1, 0)); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { for (int i = 0; i < Length; i++) { forces[d] += ResultsOfIntegration[i, 0]; } } ).Execute(); } for (int i = 0; i < D; i++) { forces[i] = MPI.Wrappers.MPIExtensions.MPISum(forces[i]); } return(forces); }
/// <summary> /// Mass defect related to divergence in continuity equation. /// </summary> /// <param name="Divergence"></param> /// <param name="Velocity"></param> /// <param name="MassDefect"></param> public static void CalculateMassDefect_Divergence(SIMPLEOperator[] Divergence, VectorField <SinglePhaseField> Velocity, SinglePhaseField MassDefect) { Debug.Assert(Divergence.Length == Velocity.Dim); for (int comp = 0; comp < Divergence.Length; comp++) { Divergence[comp].OperatorMatrix.SpMVpara(1.0, Velocity[comp].CoordinateVector, 1.0, MassDefect.CoordinateVector); MassDefect.CoordinateVector.Acc(1.0, Divergence[comp].OperatorAffine); } }
/// <summary> /// Calculates the Torque around the center of mass /// </summary> /// <param name="U"></param> /// <param name="P"></param> /// <param name="muA"></param> /// <param name="particleRadius"></param> /// <returns></returns> static public double GetTorque(VectorField <SinglePhaseField> U, SinglePhaseField P, LevelSetTracker LsTrk, double muA, double particleRadius) { var _LsTrk = LsTrk; int D = _LsTrk.GridDat.SpatialDimension; var UA = U.ToArray(); //if (D > 2) throw new NotImplementedException("Currently only 2D cases supported"); int RequiredOrder = U[0].Basis.Degree * 3 + 2; //if (RequiredOrder > agg.HMForder) // throw new ArgumentException(); Console.WriteLine("Torque coeff: {0}, order = {1}", LsTrk.CutCellQuadratureType, RequiredOrder); ConventionalDGField pA = null; double force = new double(); pA = P; ScalarFunctionEx ErrFunc = delegate(int j0, int Len, NodeSet Ns, MultidimensionalArray result) { int K = result.GetLength(1); // No nof Nodes MultidimensionalArray Grad_UARes = MultidimensionalArray.Create(Len, K, D, D);; MultidimensionalArray pARes = MultidimensionalArray.Create(Len, K); // Evaluate tangential velocity to level-set surface var Normals = _LsTrk.DataHistories[0].Current.GetLevelSetNormals(Ns, j0, Len); for (int i = 0; i < D; i++) { UA[i].EvaluateGradient(j0, Len, Ns, Grad_UARes.ExtractSubArrayShallow(-1, -1, i, -1), 0, 1); } pA.Evaluate(j0, Len, Ns, pARes); for (int j = 0; j < Len; j++) { for (int k = 0; k < K; k++) { double acc = 0.0; double acc2 = 0.0; // Calculate the torque around a circular particle with a given radius (Paper Wan and Turek 2005) acc += pARes[j, k] * Normals[j, k, 0]; acc -= (2 * muA) * Grad_UARes[j, k, 0, 0] * Normals[j, k, 0]; acc -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 1]; acc -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 1]; acc *= -Normals[j, k, 1] * particleRadius; acc2 += pARes[j, k] * Normals[j, k, 1]; acc2 -= (2 * muA) * Grad_UARes[j, k, 1, 1] * Normals[j, k, 1]; acc2 -= (muA) * Grad_UARes[j, k, 1, 0] * Normals[j, k, 0]; acc2 -= (muA) * Grad_UARes[j, k, 0, 1] * Normals[j, k, 0]; acc2 *= Normals[j, k, 0] * particleRadius; result[j, k] = acc + acc2; } } }; var SchemeHelper = LsTrk.GetXDGSpaceMetrics(new[] { LsTrk.GetSpeciesId("A") }, RequiredOrder, 1).XQuadSchemeHelper; //var SchemeHelper = new XQuadSchemeHelper(_LsTrk, momentFittingVariant, _LsTrk.GetSpeciesId("A")); CellQuadratureScheme cqs = SchemeHelper.GetLevelSetquadScheme(0, _LsTrk.Regions.GetCutCellMask()); CellQuadrature.GetQuadrature(new int[] { 1 }, _LsTrk.GridDat, cqs.Compile(_LsTrk.GridDat, RequiredOrder), delegate(int i0, int Length, QuadRule QR, MultidimensionalArray EvalResult) { ErrFunc(i0, Length, QR.Nodes, EvalResult.ExtractSubArrayShallow(-1, -1, 0)); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { for (int i = 0; i < Length; i++) { force += ResultsOfIntegration[i, 0]; } } ).Execute(); return(force); }
/// <summary> /// still a hack... /// </summary> static object InstantiateFromAttribute(FieldInfo f, InstantiateFromControlFileAttribute at, Type type, ICollection <DGField> IOFields, ICollection <DGField> RegisteredFields, //AppControl ctrl, IDictionary <string, FieldOpts> FieldOptions, GridData ctx, LevelSetTracker lstrk) { // create instance // =============== object member_value = null; Type HistoryType = null; Type VectorType = null; Type ComponentType = null; GetTypes(f, out HistoryType, out VectorType, out ComponentType); if (VectorType != null) { // vector field branch // +++++++++++++++++++ if (at.m_IsScalarField) { throw new ApplicationException("illegal use of 'InstantiateFromControlFileAttribute' (with Scalar Declaration) on a Vector class"); } int D = ctx.SpatialDimension; string[] cName = at.GetControlFileNames(f, D); string[] iName = at.GetInCodeIdentifications(f, D); // determine DG polynomial degree of basis int[] Deg = new int[D]; for (int d = 0; d < D; d++) { Deg[d] = GetDegree(cName[d], iName[d], FieldOptions); } if (at.m_DegreesMustBeEqual) { int deg0 = Deg[0]; for (int d = 1; d < D; d++) { if (Deg[d] != deg0) { StringWriter errMsg = new StringWriter(); errMsg.Write("DG Polynomial degree of fields {"); for (int dd = 0; dd < D; dd++) { errMsg.Write(cName[dd]); if (dd < D - 1) { errMsg.Write(", "); } } errMsg.Write("} must be equal, but found {"); for (int dd = 0; dd < D; dd++) { errMsg.Write(Deg[dd]); if (dd < D - 1) { errMsg.Write(", "); } } errMsg.Write("} in control file."); throw new ApplicationException(errMsg.ToString()); } } } // create instance: components SinglePhaseField[] fld = new SinglePhaseField[D]; XDGField[] xfld = new XDGField[D]; DGField[] _fld = new DGField[D]; for (int d = 0; d < D; d++) { if (ComponentType == typeof(SinglePhaseField)) { fld[d] = new SinglePhaseField(new Basis(ctx, Deg[d]), iName[d]); _fld[d] = fld[d]; } else if (ComponentType == typeof(XDGField)) { xfld[d] = new XDGField(new XDGBasis(lstrk, Deg[d]), iName[d]); _fld[d] = xfld[d]; fld = null; } else { throw new NotSupportedException("unknown type."); } RegisteredFields.Add(_fld[d]); } // create instance: Vector-Field container var ci = VectorType.GetConstructor(new Type[] { ComponentType.MakeArrayType() }); member_value = ci.Invoke(new object[] { (fld != null) ? ((object)fld) : ((object)xfld) }); //member_value = ci.Invoke( new object[] { ((object)fld) }); // io for (int d = 0; d < D; d++) { AddToIO(iName[d], _fld[d], FieldOptions, IOFields, at); } } else { // scalar field branch // +++++++++++++++++++ if (at.m_IsVectorField) { throw new ApplicationException("illegal use of 'InstantiateFromControlFileAttribute' (with Vector Declaration) on a Non-Vector class"); } // identification string cName = at.GetControlFileName(f); string iName = at.GetInCodeIdentification(f); // create basis int Deg = GetDegree(cName, iName, FieldOptions); Basis b = new Basis(ctx, Deg); // create instance DGField fld = null; if (ComponentType == typeof(SinglePhaseField)) { fld = new SinglePhaseField(b, iName); } else if (ComponentType == typeof(LevelSet)) { fld = new LevelSet(b, iName); } else if (ComponentType == typeof(XDGField)) { fld = new XDGField(new XDGBasis(lstrk, Deg), iName); } else { throw new NotImplementedException(); } RegisteredFields.Add(fld); AddToIO(iName, fld, FieldOptions, IOFields, at); member_value = fld; } // History, if desired if (HistoryType != null) { member_value = (HistoryType.GetConstructors()[0]).Invoke(new object[] { member_value }); } return(member_value); }
/// <summary> /// L2 error of some quantity derived from the state vector (e.g., /// entropy) with respect to given reference solution. The quadrature /// is determined from the settings in <see cref="IBMControl"/> /// </summary> /// <param name="quantityOfInterest"></param> /// <param name="referenceSolution"></param> /// <returns></returns> public static Query L2Error(Func <StateVector, double> quantityOfInterest, Func <double[], double, double> referenceSolution) { return(delegate(IApplication <AppControl> app, double time) { IProgram <CNSControl> program = app as IProgram <CNSControl>; if (program == null) { throw new Exception(); } ImmersedSpeciesMap speciesMap = program.SpeciesMap as ImmersedSpeciesMap; IBMControl control = program.Control as IBMControl; if (speciesMap == null || control == null) { throw new Exception( "Query is only valid for immersed boundary runs"); } SpeciesId species = speciesMap.Tracker.GetSpeciesId(control.FluidSpeciesName); int order = control.LevelSetQuadratureOrder; CellQuadratureScheme scheme = speciesMap.QuadSchemeHelper.GetVolumeQuadScheme( species, true, speciesMap.SubGrid.VolumeMask); var composititeRule = scheme.Compile(program.GridData, order); IChunkRulePair <QuadRule>[] chunkRulePairs = composititeRule.ToArray(); DGField density = program.WorkingSet.Density; VectorField <DGField> momentum = program.WorkingSet.Momentum; DGField energy = program.WorkingSet.Energy; // Construct dummy field since L2Error is currently only supported // for Field's; However, _avoid_ a projection. DGField dummy = new SinglePhaseField(new Basis(program.GridData, 0)); Material material = speciesMap.GetMaterial(double.NaN); int index = 0; double value = dummy.LxError( (ScalarFunctionEx) delegate(int j0, int Len, NodeSet nodes, MultidimensionalArray result) { MultidimensionalArray input = program.GridData.GlobalNodes.GetValue_Cell(nodes, j0, Len); Chunk chunk = chunkRulePairs[index].Chunk; QuadRule rule = chunkRulePairs[index].Rule; if (chunk.i0 != j0 || chunk.Len != Len) { throw new Exception(); } if (rule.NoOfNodes != nodes.GetLength(0)) { throw new Exception(); } MultidimensionalArray rho = MultidimensionalArray.Create(chunk.Len, rule.NoOfNodes); density.Evaluate(chunk.i0, chunk.Len, nodes, rho); MultidimensionalArray[] m = new MultidimensionalArray[CNSEnvironment.NumberOfDimensions]; for (int d = 0; d < CNSEnvironment.NumberOfDimensions; d++) { m[d] = MultidimensionalArray.Create(chunk.Len, rule.NoOfNodes); momentum[d].Evaluate(chunk.i0, chunk.Len, nodes, m[d]); } MultidimensionalArray rhoE = MultidimensionalArray.Create(chunk.Len, rule.NoOfNodes); energy.Evaluate(chunk.i0, chunk.Len, nodes, rhoE); double[] X = new double[CNSEnvironment.NumberOfDimensions]; Vector3D mVec = new Vector3D(); for (int i = 0; i < chunk.Len; i++) { for (int j = 0; j < rule.NoOfNodes; j++) { for (int d = 0; d < CNSEnvironment.NumberOfDimensions; d++) { X[d] = input[i, j, d]; mVec[d] = m[d][i, j]; } StateVector state = new StateVector(material, rho[i, j], mVec, rhoE[i, j]); double qoi = quantityOfInterest(state); Debug.Assert( !double.IsNaN(qoi), "Encountered node with unphysical state" + " (not able to determine quantity of interest)"); result[i, j] = qoi - referenceSolution(X, time); } } index++; }, (X, a, b) => (a - b) * (a - b), composititeRule); // No value is NaN, but the results. How can this be? // => All values around 0, but values in void region are a little // farther away from the exact solution // => weights in the void zone sum up to something slightly negative Debug.Assert( value >= 0, "Encountered unphysical norm even though individual values where valid." + " This indicates a problem with cut-cell quadrature."); return Math.Sqrt(value); }); }