void TerminationCriterion::accumulate_inner( PatchData& pd, OFEvaluator& of_eval, MsqError& err ) { double of_value = 0; if (terminationCriterionFlag & GRAD_FLAGS) { mGrad.resize( pd.num_free_vertices() ); bool b = of_eval.evaluate(pd, of_value, mGrad, err); MSQ_ERRRTN(err); if (!b) { MSQ_SETERR(err)("Initial patch is invalid for gradient compuation.", MsqError::INVALID_MESH); return; } } else if (terminationCriterionFlag & OF_FLAGS) { bool b = of_eval.evaluate(pd, of_value, err); MSQ_ERRRTN(err); if (!b) { MSQ_SETERR(err)("Invalid patch passed to TerminationCriterion.", MsqError::INVALID_MESH); return; } } accumulate_inner( pd, of_value, mGrad.empty() ? 0 : arrptr(mGrad), err ); MSQ_CHKERR(err); }
void TerminationCriterionTest::test_gradient_common( bool absolute, bool L2 ) { MsqPrintError err(std::cout); PatchData pd; create_twelve_hex_patch( pd, err ); ASSERT_NO_ERROR(err); const double LIMIT = 1e-4; TerminationCriterion tc; if (absolute) { if (L2) tc.add_absolute_gradient_L2_norm( LIMIT ); else tc.add_absolute_gradient_inf_norm( LIMIT ); } else { if (L2) tc.add_relative_gradient_L2_norm( LIMIT ); else tc.add_relative_gradient_inf_norm( LIMIT ); } double (*func_ptr)(const Vector3D*, int) = L2 ? &lenfunc : &maxfunc; double junk, value = 1; objFunc.mGrad = Vector3D(value,value,value); tc.reset_inner( pd, ofEval, err ); ASSERT_NO_ERROR(err); tc.reset_patch( pd, err ); ASSERT_NO_ERROR(err); std::vector<Vector3D> grad; ofEval.evaluate( pd, junk, grad, err ); ASSERT_NO_ERROR(err); double limit = LIMIT; if (!absolute) limit *= func_ptr(arrptr(grad),pd.num_free_vertices()); while (func_ptr(arrptr(grad),pd.num_free_vertices()) > limit) { CPPUNIT_ASSERT(!tc.terminate()); value *= 0.1; objFunc.mGrad = Vector3D(value,value,value); ofEval.evaluate( pd, junk, grad, err ); ASSERT_NO_ERROR(err); tc.accumulate_inner( pd, 0.0, arrptr(grad), err ); ASSERT_NO_ERROR(err); tc.accumulate_patch( pd, err ); ASSERT_NO_ERROR(err); } CPPUNIT_ASSERT(tc.terminate()); }
/*!This function checks the culling method criterion supplied to the object by the user. If the user does not supply a culling method criterion, the default criterion is NONE, and in that case, no culling is performed. If the culling method criterion is satisfied, the interior vertices of the given patch are flagged as soft_fixed. Otherwise, the soft_fixed flag is removed from each of the vertices in the patch (interior and boundary vertices). Also, if the criterion was satisfied, then the function returns true. Otherwise, the function returns false. */ bool TerminationCriterion::cull_vertices(PatchData &pd, OFEvaluator& of_eval, MsqError &err) { //PRINT_INFO("CULLING_METHOD FLAG = %i",cullingMethodFlag); //cull_bool will be changed to true if the criterion is satisfied bool b, cull_bool=false; double prev_m, init_m; switch(cullingMethodFlag){ //if no culling is requested, always return false case NONE: return cull_bool; //if culling on quality improvement absolute case QUALITY_IMPROVEMENT_ABSOLUTE: //get objective function value b = of_eval.evaluate(pd, currentOFValue, err); if (MSQ_CHKERR(err)) return false; if (!b) { MSQ_SETERR(err)(MsqError::INVALID_MESH); return false; } //if the improvement was enough, cull if(currentOFValue <= cullingEps) { cull_bool=true; } //PRINT_INFO("\ncurrentOFValue = %f, bool = %i\n",currentOFValue,cull_bool); break; //if culing on quality improvement relative case QUALITY_IMPROVEMENT_RELATIVE: //get objective function value b = of_eval.evaluate(pd, currentOFValue, err); if (MSQ_CHKERR(err)) return false; if(!b){ MSQ_SETERR(err)(MsqError::INVALID_MESH); return false; } //if the improvement was enough, cull if((currentOFValue-lowerOFBound)<= (cullingEps*(initialOFValue-lowerOFBound))) { cull_bool=true; } break; //if culling on vertex movement absolute case VERTEX_MOVEMENT_ABSOLUTE: case VERTEX_MOVEMENT_ABS_EDGE_LENGTH: //if movement was enough, cull prev_m = pd.get_max_vertex_movement_squared(previousVerticesMemento,err); MSQ_ERRZERO(err); if(prev_m <= cullingEps*cullingEps){ cull_bool=true; } break; //if culling on vertex movement relative case VERTEX_MOVEMENT_RELATIVE: //if movement was small enough, cull prev_m = pd.get_max_vertex_movement_squared(previousVerticesMemento,err); MSQ_ERRZERO(err); init_m = pd.get_max_vertex_movement_squared(initialVerticesMemento,err); MSQ_ERRZERO(err); if(prev_m <= (cullingEps*cullingEps * init_m)){ cull_bool=true; } break; case UNTANGLED_MESH: if (!patchInvertedCount) cull_bool = true; break; default: MSQ_SETERR(err)("Requested culling method not yet implemented.", MsqError::NOT_IMPLEMENTED); return false; }; //Now actually have patch data cull vertices if(cull_bool) { pd.set_free_vertices_soft_fixed(err); MSQ_ERRZERO(err); } else { pd.set_all_vertices_soft_free(err); MSQ_ERRZERO(err); } return cull_bool; }
/*!Reset function using using a PatchData object. This function is called for the inner-stopping criterion directly from the loop over mesh function in VertexMover. For outer criterion, it is called from the reset function which takes a MeshSet object. This function prepares the object to be used by setting the initial values of some of the data members. As examples, if needed, it resets the cpu timer to zero, the iteration counter to zero, and the initial and previous objective function values to the current objective function value for this patch. The return value for this function is similar to that of terminate(). The function returns false if the checked criteria have not been satisfied, and true if they have been. reset() only checks the GRADIENT_INF_NORM_ABSOLUTE, GRADIENT_L2_NORM_ABSOLUTE, and the QUALITY_IMPROVEMENT_ABSOLUTE criteria. Checking these criteria allows the QualityImprover to skip the entire optimization if the initial mesh satisfies the appropriate conditions. */ void TerminationCriterion::reset_inner(PatchData &pd, OFEvaluator& obj_eval, MsqError &err) { const unsigned long totalFlag = terminationCriterionFlag | cullingMethodFlag; // clear flag for BOUNDED_VERTEX_MOVEMENT vertexMovementExceedsBound = 0; // Use -1 to denote that this isn't initialized yet. // As all valid values must be >= 0.0, a negative // value indicates that it is uninitialized and is // always less than any valid value. maxSquaredMovement = -1; // Clear the iteration count. iterationCounter = 0; //reset the inner timer if needed if(totalFlag & CPU_TIME){ mTimer.reset(); } //GRADIENT currentGradInfNorm = initialGradInfNorm = 0.0; currentGradL2NormSquared = initialGradL2NormSquared = 0.0; if(totalFlag & GRAD_FLAGS) { if (!obj_eval.have_objective_function()) { MSQ_SETERR(err)("Error termination criteria set which uses objective " "functions, but no objective function is available.", MsqError::INVALID_STATE); return; } int num_vertices=pd.num_free_vertices(); mGrad.resize( num_vertices ); //get gradient and make sure it is valid bool b = obj_eval.evaluate(pd, currentOFValue, mGrad, err); MSQ_ERRRTN(err); if (!b) { MSQ_SETERR(err)("Initial patch is invalid for gradient computation.", MsqError::INVALID_STATE); return; } //get the gradient norms if (totalFlag & (GRADIENT_INF_NORM_ABSOLUTE|GRADIENT_INF_NORM_RELATIVE)) { currentGradInfNorm = initialGradInfNorm = Linf(mGrad); MSQ_DBGOUT_P0_ONLY(debugLevel) << par_string() << " o Initial gradient Inf norm: " << " " << RPM(initialGradInfNorm) << std::endl; } if (totalFlag & (GRADIENT_L2_NORM_ABSOLUTE|GRADIENT_L2_NORM_RELATIVE)) { currentGradL2NormSquared = initialGradL2NormSquared = length_squared(mGrad); MSQ_DBGOUT_P0_ONLY(debugLevel) << par_string() << " o Initial gradient L2 norm: " << " " << RPM(std::sqrt(initialGradL2NormSquared)) << std::endl; } //the OFvalue comes for free, so save it previousOFValue=currentOFValue; initialOFValue=currentOFValue; } //find the initial objective function value if needed and not already //computed. If we needed the gradient, we have the OF value for free. // Also, if possible, get initial OF value if writing plot file. Solvers // often supply the OF value for subsequent iterations so by calculating // the initial value we can generate OF value plots. else if ((totalFlag & OF_FLAGS) || (plotFile.is_open() && pd.num_free_vertices() && obj_eval.have_objective_function())) { //ensure the obj_ptr is not null if(!obj_eval.have_objective_function()){ MSQ_SETERR(err)("Error termination criteria set which uses objective " "functions, but no objective function is available.", MsqError::INVALID_STATE); return; } bool b = obj_eval.evaluate(pd, currentOFValue, err); MSQ_ERRRTN(err); if (!b){ MSQ_SETERR(err)("Initial patch is invalid for evaluation.",MsqError::INVALID_STATE); return; } //std::cout<<"\nReseting initial of value = "<<initialOFValue; previousOFValue=currentOFValue; initialOFValue=currentOFValue; } if (totalFlag & (GRAD_FLAGS|OF_FLAGS)) MSQ_DBGOUT_P0_ONLY(debugLevel) << par_string() << " o Initial OF value: " << " " << RPM(initialOFValue) << std::endl; // Store current vertex locations now, because we'll // need them later to compare the current movement with. if (totalFlag & VERTEX_MOVEMENT_RELATIVE) { if (initialVerticesMemento) { pd.recreate_vertices_memento( initialVerticesMemento, err ); } else { initialVerticesMemento = pd.create_vertices_memento( err ); } MSQ_ERRRTN(err); maxSquaredInitialMovement = DBL_MAX; } else { maxSquaredInitialMovement = 0; } if (terminationCriterionFlag & UNTANGLED_MESH) { globalInvertedCount = count_inverted( pd, err ); //if (innerOuterType==TYPE_OUTER) MSQ_DBGOUT_P0_ONLY(debugLevel) << par_string() << " o Num Inverted: " << " " << globalInvertedCount << std::endl; patchInvertedCount = 0; MSQ_ERRRTN(err); } if (timeStepFileType) { // If didn't already calculate gradient abive, calculate it now. if (!(totalFlag & GRAD_FLAGS)) { mGrad.resize( pd.num_free_vertices() ); obj_eval.evaluate(pd, currentOFValue, mGrad, err); err.clear(); } write_timestep( pd, mGrad.empty() ? 0 : arrptr(mGrad), err); } if (plotFile.is_open()) { // two newlines so GNU plot knows that we are starting a new data set plotFile << std::endl << std::endl; // write column headings as comment in data file plotFile << "#Iter\tCPU\tObjFunc\tGradL2\tGradInf\tMovement\tInverted" << std::endl; // write initial values plotFile << 0 << '\t' << mTimer.since_birth() << '\t' << initialOFValue << '\t' << std::sqrt( currentGradL2NormSquared ) << '\t' << currentGradInfNorm << '\t' << 0.0 << '\t' << globalInvertedCount << std::endl; } }