void ComputeL2Error(double PhysTime) { Console.WriteLine("Phystime = " + PhysTime); if ((this.Control.CircleRadius != null) != (this.Control.CutCellQuadratureType == XQuadFactoryHelper.MomentFittingVariants.ExactCircle)) { throw new ApplicationException("Illegal HMF configuration."); } if (this.Control.CircleRadius != null) { ExactCircleLevelSetIntegration.RADIUS = new double[] { this.Control.CircleRadius(PhysTime) }; } int order = Math.Max(this.u.Basis.Degree * 3, 3); XQuadSchemeHelper schH = LsTrk.GetXDGSpaceMetrics(this.LsTrk.SpeciesIdS.ToArray(), order).XQuadSchemeHelper; var uNum_A = this.u.GetSpeciesShadowField("A"); var uNum_B = this.u.GetSpeciesShadowField("B"); double uA_Err = uNum_A.L2Error(this.Control.uA_Ex.Vectorize(PhysTime), order, schH.GetVolumeQuadScheme(this.LsTrk.GetSpeciesId("A"))); double uB_Err = uNum_B.L2Error(this.Control.uB_Ex.Vectorize(PhysTime), order, schH.GetVolumeQuadScheme(this.LsTrk.GetSpeciesId("B"))); Func <double[], double, double> uJmp_Ex = ((X, t) => this.Control.uA_Ex(X, t) - this.Control.uB_Ex(X, t)); SinglePhaseField uNumJump = new SinglePhaseField(uNum_A.Basis, "Jump"); var CC = LsTrk.Regions.GetCutCellMask(); uNumJump.Acc(+1.0, uNum_A, CC); uNumJump.Acc(-1.0, uNum_B, CC); double JmpL2Err = uNumJump.L2Error(uJmp_Ex.Vectorize(PhysTime), order, schH.GetLevelSetquadScheme(0, CC)); base.QueryHandler.ValueQuery("uA_Err", uA_Err); base.QueryHandler.ValueQuery("uB_Err", uB_Err); base.QueryHandler.ValueQuery("uJmp_Err", JmpL2Err); Console.WriteLine("L2-err at t = {0}, bulk: {1}", PhysTime, Math.Sqrt(uA_Err.Pow2() + uB_Err.Pow2())); Console.WriteLine("L2-err at t = {0}, species A: {1}", PhysTime, uA_Err); Console.WriteLine("L2-err at t = {0}, species B: {1}", PhysTime, uB_Err); Console.WriteLine("L2-err at t = {0}, Jump: {1}", PhysTime, JmpL2Err); double uA_min, uA_max, uB_min, uB_max; int dummy1, dummy2; uNum_A.GetExtremalValues(out uA_min, out uA_max, out dummy1, out dummy2, this.LsTrk.Regions.GetSpeciesMask("A")); uNum_B.GetExtremalValues(out uB_min, out uB_max, out dummy1, out dummy2, this.LsTrk.Regions.GetSpeciesMask("B")); base.QueryHandler.ValueQuery("uA_Min", uA_min); base.QueryHandler.ValueQuery("uA_Max", uA_max); base.QueryHandler.ValueQuery("uB_Min", uB_min); base.QueryHandler.ValueQuery("uB_Max", uB_max); }
/// <summary> /// Determines the total curvature which is twice the mean curvature /// Please pay attention to the sign of this expression! /// </summary> /// <param name="Output">The total curvature</param> /// <param name="optionalSubGrid"> /// Subgrid which can be defined, for example for carrying out /// computations on a narrow band /// </param> /// <param name="bndMode"> /// Definition of the behavior at subgrid boundaries /// </param> public void ComputeTotalCurvatureByFlux(SinglePhaseField Output, SubGrid optionalSubGrid = null, SpatialOperator.SubGridBoundaryModes bndMode = SpatialOperator.SubGridBoundaryModes.OpenBoundary) { if (this.m_Basis.Degree <= 1) { throw new ArgumentException("For correct computation of these level set quantities, the level set has to be at least of degree 2!"); } Basis basisForNormalVec = new Basis(this.GridDat, this.m_Basis.Degree - 1); Basis basisForCurvature = new Basis(this.GridDat, this.m_Basis.Degree - 2); Func <Basis, string, SinglePhaseField> fac = (Basis b, string id) => new SinglePhaseField(b, id); VectorField <SinglePhaseField> normalVector = new VectorField <SinglePhaseField>(this.GridDat.SpatialDimension, basisForNormalVec, fac); ComputeNormalByFlux(normalVector, optionalSubGrid, bndMode); VectorField <SinglePhaseField> secondDerivatives = new VectorField <SinglePhaseField>(this.GridDat.SpatialDimension, basisForCurvature, fac); Output.Clear(); for (int i = 0; i < normalVector.Dim; i++) { secondDerivatives[i].DerivativeByFlux(1.0, normalVector[i], i, optionalSubGrid, bndMode); Output.Acc(-1.0, secondDerivatives[i], optionalSubGrid.VolumeMask); } }
/// <summary> /// Difference between <paramref name="f1"/> and <paramref name="f0"/> /// in the L2Norm. /// </summary> /// <param name="f1"></param> /// <param name="f0"></param> /// <param name="ResidualKey"> /// Optional parameter for key of residual. /// If null, the name of <paramref name="f1"/> will be used as key. /// </param> public void ComputeDifference(DGField f1, DGField f0, string ResidualKey = null) { if (f1.Basis.Degree != f0.Basis.Degree) { throw new ArgumentException("Fields must have equal degrees!"); } string key = GetDictionaryKey(Residualtype.L2Norm, f1.Identification, ResidualKey); SinglePhaseField diff = new SinglePhaseField(f1.Basis); diff.Acc(1.0, f1); diff.Acc(-1.0, f0); double res = diff.L2Norm(); m_Residuals[key] = res; }
/// <summary> /// [Multiphase] Summand for previous time steps in level-set equation. /// </summary> /// <param name="dt"></param> /// <param name="BDFOrder"></param> /// <param name="Scalar"></param> /// <param name="RhsSummand">Accumulator for the result.</param> public void ComputeRhsSummand(double dt, int BDFOrder, ScalarFieldHistory <SinglePhaseField> Scalar, SinglePhaseField RhsSummand) { for (int alpha = 1; alpha <= BDFOrder; alpha++) { RhsSummand.Acc(beta[BDFOrder - 1][alpha], Scalar[1 - alpha]); } RhsSummand.Scale(1.0 / (gamma[BDFOrder - 1] * dt)); }
/// <summary> /// [Incompressible] Summand for previous time steps in momentum equation. /// </summary> /// <param name="dt"></param> /// <param name="BDFOrder"></param> /// <param name="Velocity"></param> /// <param name="SpatialComponent">Velocity component.</param> /// <param name="RhsSummand">Accumulator for the result.</param> public void ComputeRhsSummand(double dt, int BDFOrder, VectorFieldHistory <SinglePhaseField> Velocity, int SpatialComponent, SinglePhaseField RhsSummand) { for (int alpha = 1; alpha <= BDFOrder; alpha++) { RhsSummand.Acc(beta[BDFOrder - 1][alpha], Velocity[1 - alpha][SpatialComponent]); } RhsSummand.Scale(1.0 / (gamma[BDFOrder - 1] * dt)); }
/// <summary> /// Adds an XDG field. /// </summary> public static void AddDGField(this TestingIO t, XDGField f) { var trk = f.Basis.Tracker; foreach (string spc in trk.SpeciesNames) { var fs = f.GetSpeciesShadowField(spc); SinglePhaseField fsFatClone = new SinglePhaseField(fs.Basis, fs.Identification); CellMask msk = trk.Regions.GetSpeciesMask(spc); fsFatClone.Acc(1.0, fs, msk); t.AddDGField(fsFatClone); } }
private void Filter(SinglePhaseField FiltrdField, int NoOfSweeps, CellMask CC) { Basis patchRecoveryBasis = FiltrdField.Basis; L2PatchRecovery l2pr = new L2PatchRecovery(patchRecoveryBasis, patchRecoveryBasis, CC, true); SinglePhaseField F_org = FiltrdField.CloneAs(); for (int pass = 0; pass < NoOfSweeps; pass++) { F_org.Clear(); F_org.Acc(1.0, FiltrdField); FiltrdField.Clear(); l2pr.Perform(FiltrdField, F_org); } }
/// <summary> /// Create Fields with same basis as DG Level Set /// </summary> protected override void CreateFields() { // create fields phi0 = new SinglePhaseField(DGLevSet.Basis, "phi0"); gradPhi0 = new VectorField <SinglePhaseField>(DGLevSet.GridDat.SpatialDimension.ForLoop(d => new SinglePhaseField(DGLevSet.Basis, "dPhiDG_dx[" + d + "]"))); Curvature = new SinglePhaseField(new Basis(phi0.GridDat, 0), VariableNames.Curvature); DCurvature = new SinglePhaseField(new Basis(phi0.GridDat, 0), "D" + VariableNames.Curvature); phi = new SinglePhaseField(DGLevSet.Basis, "phi"); phi.Acc(1.0, DGLevSet); mu = new SinglePhaseField(DGLevSet.Basis, "mu"); phi_Resi = new SinglePhaseField(DGLevSet.Basis, "phi_Resi"); mu_Resi = new SinglePhaseField(DGLevSet.Basis, "mu_Resi"); curvature_Resi = new SinglePhaseField(Curvature.Basis, "curvature_Resi"); // residuals: var solFields = InstantiateSolutionFields(); CurrentStateVector = new CoordinateVector(solFields); // residuals: var resFields = InstantiateResidualFields(); CurrentResidualVector = new CoordinateVector(resFields); //// Dummy Level Set //DummyLevSet = new LevelSet(new Basis(this.GridData, 1), "Levset"); //DummyLevSet.AccConstant(-1); //this.DummyLsTrk = new LevelSetTracker((GridData)(this.GridData), XQuadFactoryHelper.MomentFittingVariants.Saye, 1, new string[] { "A", "B" }, DummyLevSet); //this.DummyLsTrk.UpdateTracker(0.0); // Actual Level Set used for correction operations CorrectionLevSet = new LevelSet(phi.Basis, "Levset"); this.CorrectionLsTrk = new LevelSetTracker((GridData)(this.GridData), XQuadFactoryHelper.MomentFittingVariants.Saye, 2, new string[] { "A", "B" }, CorrectionLevSet); CorrectionLevSet.Clear(); CorrectionLevSet.Acc(1.0, phi); this.CorrectionLsTrk.UpdateTracker(0.0); // set coefficients SetCHCoefficents(); }
/// <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(new[] { "A" }, 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> /// One Iteration of the ReInitialization /// Operators must be built first /// </summary> /// <param name="ChangeRate"> /// L2-Norm of the Change-Rate in the level set in this reinit step /// </param> /// <param name="Restriction"> /// The subgrid, on which the ReInit is performed /// </param> /// <param name="IncludingInterface"> /// !! Not yet functional !! /// True, if the subgrid contains the interface, this causes all external edges of the subgrid to be treated as boundaries /// False, for the rest of the domain, thus the flux to the adjacent cells wil be evaluated /// </param> /// <returns></returns> public void ReInitSingleStep(out double ChangeRate, SubGrid Restriction = null, bool IncludingInterface = true) { if (!IncludingInterface) { throw new NotImplementedException("Untested, not yet functional!"); } using (new FuncTrace()) { /// Init Residuals Residual.Clear(); Residual.Acc(1.0, Phi); OldPhi.Clear(); OldPhi.Acc(1.0, Phi); NewPhi.Clear(); NewPhi.Acc(1.0, Phi); CellMask RestrictionMask = Restriction == null ? null : Restriction.VolumeMask; //if (Control.Upwinding && UpdateDirection && IterationCounter % 10 == 0) { if (false && Control.Upwinding && UpdateDirection) { //if (Control.Upwinding && UpdateDirection) { UpdateBulkMatrix(Restriction); } UpdateDirection = false; // RHS part RHSField.CoordinateVector.Clear(); //Operator_RHS.Evaluate(NewPhi.Mapping, RHSField.Mapping); Operator_RHS.Evaluate(double.NaN, IncludingInterface ? Restriction : null, IncludingInterface ? SubGridBoundaryModes.BoundaryEdge : SubGridBoundaryModes.InnerEdge, ArrayTools.Cat(new DGField[] { Phi }, parameterFields, new DGField[] { RHSField })); #if DEBUG RHSField.CheckForNanOrInf(); #endif // solve // ===== double[] RHS = OpAffine.CloneAs(); RHS.ScaleV(-1.0); RHS.AccV(1.0, RHSField.CoordinateVector); SolverResult Result; if (Restriction != null) { SubRHS.Clear(); SubSolution.Clear(); SubRHS.AccV(1.0, RHS, default(int[]), SubVecIdx); SubSolution.AccV(1.0, NewPhi.CoordinateVector, default(int[]), SubVecIdx); Result = slv.Solve(SubSolution, SubRHS); NewPhi.Clear(RestrictionMask); NewPhi.CoordinateVector.AccV(1.0, SubSolution, SubVecIdx, default(int[])); } else { Result = slv.Solve(NewPhi.CoordinateVector, RHS); } #if Debug OpMatrix.SpMV(-1.0, NewPhi.CoordinateVector, 1.0, RHS); Console.WriteLine("LinearSolver: {0} Iterations, Converged={1}, Residual = {2} ", Result.NoOfIterations, Result.Converged, RHS.L2Norm()); #endif // Apply underrelaxation Phi.Clear(RestrictionMask); Phi.Acc(1 - underrelaxation, OldPhi, RestrictionMask); Phi.Acc(underrelaxation, NewPhi, RestrictionMask); Residual.Acc(-1.0, Phi, RestrictionMask); ChangeRate = Residual.L2Norm(RestrictionMask); //Calculate LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi, RestrictionMask); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient, RestrictionMask); if (Control.Upwinding) { //RestrictionMask.GetBitMask(); for (int i = 0; i < MeanLevelSetGradient.CoordinateVector.Length; i++) { NewDirection[i] = Math.Sign(MeanLevelSetGradient.CoordinateVector[i]); //NewDirection[i] = MeanLevelSetGradient.CoordinateVector[i]; OldDirection[i] -= NewDirection[i]; } double MaxDiff = OldDirection.L2Norm(); //if (MaxDiff > 1E-20 && IterationCounter % 10 == 0 ) { //if (MaxDiff > 1E-20) { // Console.WriteLine("Direction Values differ by {0}", MaxDiff); if (MaxDiff > 0.2) { //UpdateDirection = true; //Console.WriteLine("Direction Values differ by {0} => Updating ReInit-Matrix", MaxDiff); } ; //} //Console.WriteLine("HACK!!! Updating Upwind Matrix everytime!"); //UpdateDirection = true; // Reset Value OldDirection.Clear(); OldDirection.AccV(1.0, NewDirection); } } }
// Local Variables for Iteration // <summary> // Counter for Iteration Steps // </summary> //double OldResidual = double.MaxValue; //int divergencecounter = 0; ///// <summary> ///// Checks for Reaching Max. Number of Iterations and Divergence of Algorithm ///// </summary> ///// <param name="Residual">Change Rate of the Algorithm</param> ///// <returns>Reaching Max Iterations, Aborts when diverged</returns> //public bool CheckAbortCriteria(double Residual, int IterationCounter) { // if (Residual <= ConvergenceCriterion) { // Console.WriteLine("EllipticReInit converged after {0} Iterations ", IterationCounter); // return true; // } // if (Residual >= OldResidual) divergencecounter++; // else divergencecounter = 0; // if (IterationCounter >= MaxIteration) { // Console.WriteLine("Elliptic Reinit Max Iterations Reached"); // return true; // }; // if (divergencecounter > MaxIteration / 2) { // Console.WriteLine("Elliptic Reinit diverged - Aborting"); // throw new ApplicationException(); // } // OldResidual = Residual; // IterationCounter++; // return false; //} //bool PreviouslyOnSubgrid = false; /// <summary> /// Updates the Operator Matrix after level-set motion /// </summary> /// <param name="Restriction"> /// The subgrid, on which the ReInit is performed /// </param> /// <param name="IncludingInterface"> /// !! Not yet functional !! /// True, if the subgrid contains the interface, this causes all external edges of the subgrid to be treated as boundaries /// False, for the rest of the domain, thus the flux to the adjacent cells wil be evaluated /// </param> public void UpdateOperators(SubGrid Restriction = null, bool IncludingInterface = true) { if (!IncludingInterface) { throw new NotImplementedException("Untested, not yet functional!"); } using (new FuncTrace()) { //using (var slv = new ilPSP.LinSolvers.MUMPS.MUMPSSolver()) { //using (var slv = new ilPSP.LinSolvers.PARDISO.PARDISOSolver()) { //using (var slv = new ilPSP.LinSolvers.HYPRE.GMRES()) { if (Control.Upwinding) { OldPhi.Clear(); OldPhi.Acc(1.0, Phi); //Calculate LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi, Restriction?.VolumeMask); //LevelSetGradient.Gradient(1.0, Phi); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient, Restriction?.VolumeMask); //MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient); } if (slv != null) { slv.Dispose(); } slv = Control.solverFactory(); OpMatrix_interface.Clear(); OpAffine_interface.Clear(); // Build the Quadrature-Scheme for the interface operator // Note: The HMF-Quadrature over a surface is formally a volume quadrature, since it uses the volume quadrature nodes. //XSpatialOperatorExtensions.ComputeMatrixEx(Operator_interface, ////Operator_interface.ComputeMatrixEx( // LevelSetTracker, // Phi.Mapping, // null, // Phi.Mapping, // OpMatrix_interface, // OpAffine_interface, // false, // 0, // false, // subGrid:Restriction, // whichSpc: LevelSetTracker.GetSpeciesId("A") // ); XSpatialOperatorMk2.XEvaluatorLinear mtxBuilder = Operator_interface.GetMatrixBuilder(LevelSetTracker, Phi.Mapping, null, Phi.Mapping); MultiphaseCellAgglomerator dummy = LevelSetTracker.GetAgglomerator(LevelSetTracker.SpeciesIdS.ToArray(), Phi.Basis.Degree * 2 + 2, 0.0); //mtxBuilder.SpeciesOperatorCoefficients[LevelSetTracker.GetSpeciesId("A")].CellLengthScales = dummy.CellLengthScales[LevelSetTracker.GetSpeciesId("A")]; mtxBuilder.CellLengthScales.Add(LevelSetTracker.GetSpeciesId("A"), dummy.CellLengthScales[LevelSetTracker.GetSpeciesId("A")]); mtxBuilder.time = 0; mtxBuilder.MPITtransceive = false; mtxBuilder.ComputeMatrix(OpMatrix_interface, OpAffine_interface); // Regenerate OpMatrix for subgrid -> adjacent cells must be trated as boundary if (Restriction != null) { OpMatrix_bulk.Clear(); OpAffine_bulk.Clear(); //Operator_bulk.ComputeMatrix( // Phi.Mapping, // parameterFields, // Phi.Mapping, // OpMatrix_bulk, OpAffine_bulk, // OnlyAffine: false, sgrd: Restriction); EdgeQuadratureScheme edgescheme; //if (Control.Upwinding) { // edgescheme = new EdgeQuadratureScheme(true, IncludingInterface ? Restriction.AllEdgesMask : null); //} //else { edgescheme = new EdgeQuadratureScheme(true, IncludingInterface ? Restriction.InnerEdgesMask : null); //} Operator_bulk.ComputeMatrixEx(Phi.Mapping, parameterFields, Phi.Mapping, OpMatrix_bulk, OpAffine_bulk, false, 0, edgeQuadScheme: edgescheme, volQuadScheme: new CellQuadratureScheme(true, IncludingInterface ? Restriction.VolumeMask : null) ); //PreviouslyOnSubgrid = true; } // recalculate full Matrix //else if (PreviouslyOnSubgrid) { else { OpMatrix_bulk.Clear(); OpAffine_bulk.Clear(); Operator_bulk.ComputeMatrixEx(Phi.Mapping, parameterFields, Phi.Mapping, OpMatrix_bulk, OpAffine_bulk, false, 0 ); //PreviouslyOnSubgrid = false; } /// Compose the Matrix /// This is symmetric due to the symmetry of the SIP and the penalty term OpMatrix.Clear(); OpMatrix.Acc(1.0, OpMatrix_bulk); OpMatrix.Acc(1.0, OpMatrix_interface); OpMatrix.AssumeSymmetric = !Control.Upwinding; //OpMatrix.AssumeSymmetric = false; /// Compose the RHS of the above operators. (-> Boundary Conditions) /// This does NOT include the Nonlinear RHS, which will be added later OpAffine.Clear(); OpAffine.AccV(1.0, OpAffine_bulk); OpAffine.AccV(1.0, OpAffine_interface); #if Debug ilPSP.Connectors.Matlab.BatchmodeConnector matlabConnector; matlabConnector = new BatchmodeConnector(); #endif if (Restriction != null) { SubVecIdx = Phi.Mapping.GetSubvectorIndices(Restriction, true, new int[] { 0 }); int L = SubVecIdx.Length; SubMatrix = new MsrMatrix(L); SubRHS = new double[L]; SubSolution = new double[L]; OpMatrix.AccSubMatrixTo(1.0, SubMatrix, SubVecIdx, default(int[]), SubVecIdx, default(int[])); slv.DefineMatrix(SubMatrix); #if Debug Console.WriteLine("ConditionNumber of ReInit-Matrix is " + SubMatrix.condest().ToString("E")); #endif } else { slv.DefineMatrix(OpMatrix); #if Debug Console.WriteLine("ConditionNumber of ReInit-Matrix is " + OpMatrix.condest().ToString("E")); #endif } } }
/// <summary> /// Computes the new level set field at time <paramref name="Phystime"/> + <paramref name="dt"/>. /// This is a 'driver function' which provides a universal interface to the various level set evolution algorithms. /// It also acts as a callback to the time stepper (see <see cref="m_BDF_Timestepper"/> resp. <see cref="m_RK_Timestepper"/>), /// i.e. it matches the signature of /// <see cref="BoSSS.Solution.XdgTimestepping.DelUpdateLevelset"/>. /// </summary> /// <param name="Phystime"></param> /// <param name="dt"></param> /// <param name="CurrentState"> /// The current solution (velocity and pressure), since the complete solution is provided by the time stepper, /// only the velocity components(supposed to be at the beginning) are used. /// </param> /// <param name="underrelax"> /// </param> public double DelUpdateLevelSet(DGField[] CurrentState, double Phystime, double dt, double underrelax, bool incremental) { using (new FuncTrace()) { //dt *= underrelax; int D = base.Grid.SpatialDimension; int iTimestep = hack_TimestepIndex; DGField[] EvoVelocity = CurrentState.GetSubVector(0, D); // ======================================================== // Backup old level-set, in order to compute the residual // ======================================================== SinglePhaseField LsBkUp = new SinglePhaseField(this.LevSet.Basis); LsBkUp.Acc(1.0, this.LevSet); CellMask oldCC = LsTrk.Regions.GetCutCellMask(); // ==================================================== // set evolution velocity, but only on the CUT-cells // ==================================================== #region Calculate density averaged Velocity for each cell ConventionalDGField[] meanVelocity = GetMeanVelocityFromXDGField(EvoVelocity); #endregion // =================================================================== // backup interface properties (mass conservation, surface changerate) // =================================================================== #region backup interface props double oldSurfVolume = 0.0; double oldSurfLength = 0.0; double SurfChangerate = 0.0; if (this.Control.CheckInterfaceProps) { oldSurfVolume = XNSEUtils.GetSpeciesArea(this.LsTrk, LsTrk.GetSpeciesId("A")); oldSurfLength = XNSEUtils.GetInterfaceLength(this.LsTrk); SurfChangerate = EnergyUtils.GetSurfaceChangerate(this.LsTrk, meanVelocity, this.m_HMForder); } #endregion // ==================================================== // perform level-set evolution // ==================================================== #region level-set evolution // set up for Strang splitting SinglePhaseField DGLevSet_old; if (incremental) { DGLevSet_old = this.DGLevSet.Current.CloneAs(); } else { DGLevSet_old = this.DGLevSet[0].CloneAs(); } // set up for underrelaxation SinglePhaseField DGLevSet_oldIter = this.DGLevSet.Current.CloneAs(); //PlotCurrentState(hack_Phystime, new TimestepNumber(new int[] { hack_TimestepIndex, 0 }), 2); // actual evolution switch (this.Control.Option_LevelSetEvolution) { case LevelSetEvolution.None: throw new ArgumentException("illegal call"); case LevelSetEvolution.FastMarching: { NarrowMarchingBand.Evolve_Mk2( dt, this.LsTrk, DGLevSet_old, this.DGLevSet.Current, this.DGLevSetGradient, meanVelocity, this.ExtensionVelocity.Current.ToArray(), //new DGField[] { LevSetSrc }, this.m_HMForder, iTimestep); //FastMarchReinitSolver = new FastMarchReinit(DGLevSet.Current.Basis); //CellMask Accepted = LsTrk.Regions.GetCutCellMask(); //CellMask ActiveField = LsTrk.Regions.GetNearFieldMask(1); //CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); //FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); break; } case LevelSetEvolution.Fourier: { Fourier_Timestepper.moveLevelSet(dt, meanVelocity); if (incremental) { Fourier_Timestepper.updateFourierLevSet(); } Fourier_LevSet.ProjectToDGLevelSet(this.DGLevSet.Current, this.LsTrk); break; } case LevelSetEvolution.Prescribed: { this.DGLevSet.Current.Clear(); this.DGLevSet.Current.ProjectField(1.0, this.Control.Phi.Vectorize(Phystime + dt)); break; } case LevelSetEvolution.ScalarConvection: { var LSM = new LevelSetMover(EvoVelocity, this.ExtensionVelocity, this.LsTrk, XVelocityProjection.CutCellVelocityProjectiontype.L2_plain, this.DGLevSet, this.BcMap); int check1 = this.ExtensionVelocity.PushCount; int check2 = this.DGLevSet.PushCount; this.DGLevSet[1].Clear(); this.DGLevSet[1].Acc(1.0, DGLevSet_old); LSM.Advect(dt); if (check1 != this.ExtensionVelocity.PushCount) { throw new ApplicationException(); } if (check2 != this.DGLevSet.PushCount) { throw new ApplicationException(); } break; } case LevelSetEvolution.ExtensionVelocity: { DGLevSetGradient.Clear(); DGLevSetGradient.Gradient(1.0, DGLevSet.Current); ExtVelMover.Advect(dt); // Fast Marching: Specify the Domains first // Perform Fast Marching only on the Far Field if (this.Control.AdaptiveMeshRefinement) { int NoCells = ((GridData)this.GridData).Cells.Count; BitArray Refined = new BitArray(NoCells); for (int j = 0; j < NoCells; j++) { if (((GridData)this.GridData).Cells.GetCell(j).RefinementLevel > 0) { Refined[j] = true; } } CellMask Accepted = new CellMask(this.GridData, Refined); CellMask AcceptedNeigh = Accepted.AllNeighbourCells(); Accepted = Accepted.Union(AcceptedNeigh); CellMask ActiveField = Accepted.Complement(); CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); } else { CellMask Accepted = LsTrk.Regions.GetNearFieldMask(1); CellMask ActiveField = Accepted.Complement(); CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); } //SubGrid AcceptedGrid = new SubGrid(Accepted); //ReInitPDE.ReInitialize(Restriction: AcceptedGrid); //CellMask ActiveField = Accepted.Complement(); //CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); //FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); //ReInitPDE.ReInitialize(); break; } default: throw new ApplicationException(); } // performing underrelaxation if (underrelax < 1.0) { this.DGLevSet.Current.Scale(underrelax); this.DGLevSet.Current.Acc((1.0 - underrelax), DGLevSet_oldIter); } //PlotCurrentState(hack_Phystime, new TimestepNumber(new int[] { hack_TimestepIndex, 1 }), 2); #endregion // ====================== // postprocessing // ======================= if (this.Control.ReInitPeriod > 0 && hack_TimestepIndex % this.Control.ReInitPeriod == 0) { Console.WriteLine("Filtering DG-LevSet"); SinglePhaseField FiltLevSet = new SinglePhaseField(DGLevSet.Current.Basis); FiltLevSet.AccLaidBack(1.0, DGLevSet.Current); Filter(FiltLevSet, 2, oldCC); DGLevSet.Current.Clear(); DGLevSet.Current.Acc(1.0, FiltLevSet); Console.WriteLine("FastMarchReInit performing FirstOrderReInit"); FastMarchReinitSolver = new FastMarchReinit(DGLevSet.Current.Basis); CellMask Accepted = LsTrk.Regions.GetCutCellMask(); CellMask ActiveField = LsTrk.Regions.GetNearFieldMask(1); CellMask NegativeField = LsTrk.Regions.GetSpeciesMask("A"); FastMarchReinitSolver.FirstOrderReinit(DGLevSet.Current, Accepted, NegativeField, ActiveField); } #region ensure continuity // make level set continuous CellMask CC = LsTrk.Regions.GetCutCellMask4LevSet(0); CellMask Near1 = LsTrk.Regions.GetNearMask4LevSet(0, 1); CellMask PosFF = LsTrk.Regions.GetLevelSetWing(0, +1).VolumeMask; ContinuityEnforcer.MakeContinuous(this.DGLevSet.Current, this.LevSet, Near1, PosFF); if (this.Control.Option_LevelSetEvolution == LevelSetEvolution.FastMarching) { CellMask Nearband = Near1.Union(CC); this.DGLevSet.Current.Clear(Nearband); this.DGLevSet.Current.AccLaidBack(1.0, this.LevSet, Nearband); //ContinuityEnforcer.SetFarField(this.DGLevSet.Current, Near1, PosFF); } //PlotCurrentState(hack_Phystime, new TimestepNumber(new int[] { hack_TimestepIndex, 2 }), 2); #endregion for (int d = 0; d < D; d++) { this.XDGvelocity.Velocity[d].UpdateBehaviour = BehaveUnder_LevSetMoovement.AutoExtrapolate; } // =============== // tracker update // =============== this.LsTrk.UpdateTracker(Phystime + dt, incremental: true); // update near field (in case of adaptive mesh refinement) if (this.Control.AdaptiveMeshRefinement && this.Control.Option_LevelSetEvolution == LevelSetEvolution.FastMarching) { Near1 = LsTrk.Regions.GetNearMask4LevSet(0, 1); PosFF = LsTrk.Regions.GetLevelSetWing(0, +1).VolumeMask; ContinuityEnforcer.SetFarField(this.DGLevSet.Current, Near1, PosFF); ContinuityEnforcer.SetFarField(this.LevSet, Near1, PosFF); } // ================================================================== // check interface properties (mass conservation, surface changerate) // ================================================================== if (this.Control.CheckInterfaceProps) { double currentSurfVolume = XNSEUtils.GetSpeciesArea(this.LsTrk, LsTrk.GetSpeciesId("A")); double massChange = ((currentSurfVolume - oldSurfVolume) / oldSurfVolume) * 100; Console.WriteLine("Change of mass = {0}%", massChange); double currentSurfLength = XNSEUtils.GetInterfaceLength(this.LsTrk); double actualSurfChangerate = (currentSurfLength - oldSurfLength) / dt; Console.WriteLine("Interface divergence = {0}", SurfChangerate); Console.WriteLine("actual surface changerate = {0}", actualSurfChangerate); } // ================== // compute residual // ================== var newCC = LsTrk.Regions.GetCutCellMask(); LsBkUp.Acc(-1.0, this.LevSet); double LevSetResidual = LsBkUp.L2Norm(newCC.Union(oldCC)); return(LevSetResidual); } }
protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { Console.WriteLine(" Timestep # " + TimestepNo + ", phystime = " + phystime); //phystime = 1.8; LsUpdate(phystime); // operator-matrix assemblieren MsrMatrix OperatorMatrix = new MsrMatrix(u.Mapping, u.Mapping); double[] Affine = new double[OperatorMatrix.RowPartitioning.LocalLength]; // Agglomerator setup MultiphaseCellAgglomerator Agg = LsTrk.GetAgglomerator(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, QuadOrder, this.THRESHOLD); // plausibility of cell length scales if (SER_PAR_COMPARISON) { TestLengthScales(QuadOrder, TimestepNo); } Console.WriteLine("Inter-Process agglomeration? " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.InterProcessAgglomeration); if (this.THRESHOLD > 0.01) { TestAgglomeration_Extraploation(Agg); TestAgglomeration_Projection(QuadOrder, Agg); } CheckExchange(true); CheckExchange(false); // operator matrix assembly XSpatialOperatorMk2.XEvaluatorLinear mtxBuilder = Op.GetMatrixBuilder(base.LsTrk, u.Mapping, null, u.Mapping); mtxBuilder.time = 0.0; mtxBuilder.ComputeMatrix(OperatorMatrix, Affine); Agg.ManipulateMatrixAndRHS(OperatorMatrix, Affine, u.Mapping, u.Mapping); // mass matrix factory var Mfact = LsTrk.GetXDGSpaceMetrics(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, QuadOrder, 1).MassMatrixFactory;// new MassMatrixFactory(u.Basis, Agg); // Mass matrix/Inverse Mass matrix //var MassInv = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, true, LsTrk.GetSpeciesId("B")); var Mass = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, false, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(Mass, default(double[]), u.Mapping, u.Mapping); var MassInv = Mass.InvertBlocks(OnlyDiagonal: true, Subblocks: true, ignoreEmptyBlocks: true, SymmetricalInversion: false); // test that operator depends only on B-species values double DepTest = LsTrk.Regions.GetSpeciesSubGrid("B").TestMatrixDependency(OperatorMatrix, u.Mapping, u.Mapping); Console.WriteLine("Matrix dependency test: " + DepTest); Assert.LessOrEqual(DepTest, 0.0); // diagnostic output Console.WriteLine("Number of Agglomerations (all species): " + Agg.TotalNumberOfAgglomerations); Console.WriteLine("Number of Agglomerations (species 'B'): " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.SourceCells.NoOfItemsLocally.MPISum()); // operator auswerten: double[] x = new double[Affine.Length]; BLAS.daxpy(x.Length, 1.0, Affine, 1, x, 1); OperatorMatrix.SpMVpara(1.0, u.CoordinateVector, 1.0, x); MassInv.SpMV(1.0, x, 0.0, du_dx.CoordinateVector); Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).Extrapolate(du_dx.Mapping); // markieren, wo ueberhaupt A und B sind Bmarker.AccConstant(1.0, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Amarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("A").VolumeMask); if (usePhi0 && usePhi1) { Xmarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("X").VolumeMask); } // compute error ERR.Clear(); ERR.Acc(1.0, du_dx_Exact, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); ERR.Acc(-1.0, du_dx, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double L2Err = ERR.L2Norm(LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Console.WriteLine("L2 Error: " + L2Err); XERR.Clear(); XERR.GetSpeciesShadowField("B").Acc(1.0, ERR, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double xL2Err = XERR.L2Norm(); Console.WriteLine("L2 Error (in XDG space): " + xL2Err); // check error double ErrorThreshold = 1.0e-1; if (this.MomentFittingVariant == XQuadFactoryHelper.MomentFittingVariants.OneStepGaussAndStokes) { ErrorThreshold = 1.0e-6; // HMF is designed for such integrands and should perform close to machine accuracy; on general integrands, the precision is different. } bool IsPassed = ((L2Err <= ErrorThreshold || this.THRESHOLD <= ErrorThreshold) && xL2Err <= ErrorThreshold); if (IsPassed) { Console.WriteLine("Test PASSED"); } else { Console.WriteLine("Test FAILED: check errors."); //PlotCurrentState(phystime, TimestepNo, 3); } if (TimestepNo > 1) { if (this.THRESHOLD > ErrorThreshold) { // without agglomeration, the error in very tiny cut-cells may be large over the whole cell // However, the error in the XDG-space should be small under all circumstances Assert.LessOrEqual(L2Err, ErrorThreshold, "DG L2 error of computing du_dx"); } Assert.LessOrEqual(xL2Err, ErrorThreshold, "XDG L2 error of computing du_dx"); } // return/Ende base.NoOfTimesteps = 17; //base.NoOfTimesteps = 2; dt = 0.3; return(dt); }
/// <summary> /// In this method, the vector normal to the surface as well as certain derivatives /// of it that we intend to use are set up /// </summary> protected void NormalVec(double deltax) { wx = new SinglePhaseField(m_gradBasis, "wx"); wy = new SinglePhaseField(m_gradBasis, "wy"); if (m_Context.Grid.SpatialDimension == 2) { m_Gradw = new VectorField <SinglePhaseField>(wx, wy); } else if (m_Context.Grid.SpatialDimension == 3) { wz = new SinglePhaseField(m_gradBasis, "wz"); m_Gradw = new VectorField <SinglePhaseField>(wx, wy, wz); } else { throw new NotSupportedException("only spatial dimension 2 and 3 are supported."); } SinglePhaseField absval = new SinglePhaseField(m_gradBasis); m_Gradw[0].Derivative(1.0, m_Field, 0); m_Gradw[1].Derivative(1.0, m_Field, 1); if (m_Context.Grid.SpatialDimension == 3) { m_Gradw[2].Derivative(1.0, m_Field, 2); } /* * Normalization of the gradient vector * Might be better implemented using a * Function */ absval.ProjectAbs(1.0, m_Gradw); m_Gradw[0].ProjectQuotient(1.0, m_Gradw[0], absval, null, false); m_Gradw[1].ProjectQuotient(1.0, m_Gradw[1], absval, null, false); if (m_Context.Grid.SpatialDimension == 3) { m_Gradw[2].ProjectQuotient(1.0, m_Gradw[2], absval, null, false); } SinglePhaseField n1x = new SinglePhaseField(m_derivsBasis, "n1x"); SinglePhaseField n2x = new SinglePhaseField(m_derivsBasis, "n2x"); SinglePhaseField n1y = new SinglePhaseField(m_derivsBasis, "n1y"); SinglePhaseField n2y = new SinglePhaseField(m_derivsBasis, "n2y"); if (m_Context.Grid.SpatialDimension == 3) { SinglePhaseField n3x = new SinglePhaseField(m_derivsBasis, "n3x"); SinglePhaseField n3y = new SinglePhaseField(m_derivsBasis, "n3y"); SinglePhaseField n1z = new SinglePhaseField(m_derivsBasis, "n1z"); SinglePhaseField n2z = new SinglePhaseField(m_derivsBasis, "n2z"); SinglePhaseField n3z = new SinglePhaseField(m_derivsBasis, "n3z"); m_normderivs = new VectorField <SinglePhaseField>(n1x, n2x, n3x, n1y, n2y, n3y, n1z, n2z, n3z); /* Partial derivatives of the normal vector. * These quantities are used in the product rule applied to the surface projection * and to compute the mean curvature */ m_normderivs[0].Derivative(1.0, m_Gradw[0], 0); m_normderivs[1].Derivative(1.0, m_Gradw[1], 0); m_normderivs[2].Derivative(1.0, m_Gradw[2], 0); m_normderivs[3].Derivative(1.0, m_Gradw[0], 1); m_normderivs[4].Derivative(1.0, m_Gradw[1], 1); m_normderivs[5].Derivative(1.0, m_Gradw[2], 1); m_normderivs[6].Derivative(1.0, m_Gradw[0], 2); m_normderivs[7].Derivative(1.0, m_Gradw[1], 2); m_normderivs[8].Derivative(1.0, m_Gradw[2], 2); } else { m_normderivs = new VectorField <SinglePhaseField>(n1x, n2x, n1y, n2y); /* Partial derivatives of the normal vector. * These quantities are used in the product rule applied to the surface projection * and to compute the mean curvature */ m_normderivs[0].Derivative(1.0, m_Gradw[0], 0); m_normderivs[1].Derivative(1.0, m_Gradw[1], 0); m_normderivs[2].Derivative(1.0, m_Gradw[0], 1); m_normderivs[3].Derivative(1.0, m_Gradw[1], 1); } /* * Divergence of the surface normal vector yields the total curvature * which is twice the mean curvature * ATTENTION: here with inverse sign!!!!! */ m_Curvature = new SinglePhaseField(m_derivsBasis); m_Curvature.Acc(1.0, m_normderivs[0]); m_Curvature.Acc(1.0, m_normderivs[4]); if (m_Context.Grid.SpatialDimension == 3) { m_Curvature.Acc(1.0, m_normderivs[8]); } m_sign = new SinglePhaseField(m_Basis, "sign"); SinglePhaseField quot = new SinglePhaseField(m_Basis, "quot"); SinglePhaseField sqrt = new SinglePhaseField(m_Basis, "sqrt"); sqrt.ProjectProduct(1.0, m_Field, m_Field); sqrt.AccConstant(0.01); quot.ProjectPow(1.0, sqrt, deltax * deltax); m_sign.ProjectQuotient(1.0, m_Field, quot); }
public void Advect(double dt) { VelocityExtender.ConstructExtension(nearfield: nearfield); LevelSet.Acc(-dt, Extension); }
protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { Console.WriteLine(" Timestep # " + TimestepNo + ", phystime = " + phystime); //phystime = 1.8; LsUpdate(phystime); // operator-matrix assemblieren MsrMatrix OperatorMatrix = new MsrMatrix(u.Mapping, u.Mapping); double[] Affine = new double[OperatorMatrix.RowPartitioning.LocalLength]; MultiphaseCellAgglomerator Agg; MassMatrixFactory Mfact; // Agglomerator setup int quadOrder = Op.QuadOrderFunction(new int[] { u.Basis.Degree }, new int[0], new int[] { u.Basis.Degree }); //Agg = new MultiphaseCellAgglomerator(new CutCellMetrics(MomentFittingVariant, quadOrder, LsTrk, ), this.THRESHOLD, false); Agg = LsTrk.GetAgglomerator(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, this.THRESHOLD); Console.WriteLine("Inter-Process agglomeration? " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.InterProcessAgglomeration); if (this.THRESHOLD > 0.01) { TestAgglomeration_Extraploation(Agg); TestAgglomeration_Projection(quadOrder, Agg); } // operator matrix assembly Op.ComputeMatrixEx(LsTrk, u.Mapping, null, u.Mapping, OperatorMatrix, Affine, false, 0.0, true, Agg.CellLengthScales, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(OperatorMatrix, Affine, u.Mapping, u.Mapping); // mass matrix factory Mfact = LsTrk.GetXDGSpaceMetrics(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, 1).MassMatrixFactory;// new MassMatrixFactory(u.Basis, Agg); // Mass matrix/Inverse Mass matrix //var MassInv = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, true, LsTrk.GetSpeciesId("B")); var Mass = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, false, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(Mass, default(double[]), u.Mapping, u.Mapping); var MassInv = Mass.InvertBlocks(OnlyDiagonal: true, Subblocks: true, ignoreEmptyBlocks: true, SymmetricalInversion: false); // test that operator depends only on B-species values double DepTest = LsTrk.Regions.GetSpeciesSubGrid("B").TestMatrixDependency(OperatorMatrix, u.Mapping, u.Mapping); Console.WriteLine("Matrix dependency test: " + DepTest); Assert.LessOrEqual(DepTest, 0.0); // diagnostic output Console.WriteLine("Number of Agglomerations (all species): " + Agg.TotalNumberOfAgglomerations); Console.WriteLine("Number of Agglomerations (species 'B'): " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.SourceCells.NoOfItemsLocally.MPISum()); // operator auswerten: double[] x = new double[Affine.Length]; BLAS.daxpy(x.Length, 1.0, Affine, 1, x, 1); OperatorMatrix.SpMVpara(1.0, u.CoordinateVector, 1.0, x); MassInv.SpMV(1.0, x, 0.0, du_dx.CoordinateVector); Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).Extrapolate(du_dx.Mapping); // markieren, wo ueberhaupt A und B sind Bmarker.AccConstant(1.0, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Amarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("A").VolumeMask); Xmarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("X").VolumeMask); // compute error ERR.Clear(); ERR.Acc(1.0, du_dx_Exact, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); ERR.Acc(-1.0, du_dx, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double L2Err = ERR.L2Norm(LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Console.WriteLine("L2 Error: " + L2Err); XERR.Clear(); XERR.GetSpeciesShadowField("B").Acc(1.0, ERR, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double xL2Err = XERR.L2Norm(); Console.WriteLine("L2 Error (in XDG space): " + xL2Err); // check error if (this.THRESHOLD > 0.01) { // without agglomeration, the error in very tiny cut-cells may be large over the whole cell // However, the error in the XDG-space should be small under all circumstances Assert.LessOrEqual(L2Err, 1.0e-6); } Assert.LessOrEqual(xL2Err, 1.0e-6); bool IsPassed = ((L2Err <= 1.0e-6 || this.THRESHOLD <= 0.01) && xL2Err <= 1.0e-7); if (IsPassed) { Console.WriteLine("Test PASSED"); } else { Console.WriteLine("Test FAILED: check errors."); } // return/Ende base.NoOfTimesteps = 17; //base.NoOfTimesteps = 2; dt = 0.3; return(dt); }
/// <summary> /// Update scalar field variables after solving scalar equation. /// </summary> /// <param name="SolverConf"></param> /// <param name="ModeRelaxScalar"></param> /// <param name="relax_scalar"></param> /// <param name="Scalar"></param> /// <param name="ScalarRes"></param> /// <param name="ScalarMean"></param> /// <param name="Rho"></param> /// <param name="Eta"></param> /// <param name="RhoMatrix"></param> /// <param name="EoS"></param> /// <param name="ThermodynamicPressure">Null for multiphase flows.</param> public static void UpdateScalarFieldVariables(SIMPLEControl SolverConf, RelaxationTypes ModeRelaxScalar, double relax_scalar, ScalarFieldHistory <SinglePhaseField> Scalar, SinglePhaseField ScalarRes, SinglePhaseField ScalarMean, SinglePhaseField Rho, SinglePhaseField Eta, QuadratureMatrix RhoMatrix, MaterialLaw EoS, SinglePhaseField ThermodynamicPressure, bool UpdateRhoVisc = true) { using (new FuncTrace()) { // Explicit Under-Relaxation of scalar variable // ============================================ if (ModeRelaxScalar == RelaxationTypes.Explicit) { // phi = alpha * phi_new + (1-alpha) * phi_old Scalar.Current.Scale(relax_scalar); Scalar.Current.Acc((1.0 - relax_scalar), ScalarRes); } // Scalar residual // =============== ScalarRes.Scale(-1.0); ScalarRes.Acc(1.0, Scalar.Current); // ScalarMean // ========== ScalarMean.Clear(); ScalarMean.AccLaidBack(1.0, Scalar.Current); // Thermodynamic pressure - only for Low-Mach number flows // ======================================================= switch (SolverConf.PhysicsMode) { case PhysicsMode.LowMach: LowMachSIMPLEControl lowMachConf = SolverConf as LowMachSIMPLEControl; if (lowMachConf.ThermodynamicPressureMode == ThermodynamicPressureMode.MassDetermined) { ThermodynamicPressure.Clear(); ThermodynamicPressure.AccConstant(((MaterialLawLowMach)EoS).GetMassDeterminedThermodynamicPressure(lowMachConf.InitialMass.Value, Scalar.Current)); } break; case PhysicsMode.Multiphase: break; default: throw new ApplicationException(); } if (UpdateRhoVisc) { // Density // ======= Rho.Clear(); Rho.ProjectFunction(1.0, EoS.GetDensity, null, Scalar.Current); RhoMatrix.Update(); // Viscosity // ========= Eta.Clear(); Eta.ProjectFunction(1.0, EoS.GetViscosity, null, Scalar.Current); } } }