00001
00002
00003
00004
00005
00006
00007 #include <fstream>
00008 #include "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010 #include "Teuchos_CommHelpers.hpp"
00011 #include "Phalanx.hpp"
00012
00013 template<typename EvalT, typename Traits>
00014 QCAD::ResponseFieldAverage<EvalT, Traits>::
00015 ResponseFieldAverage(Teuchos::ParameterList& p,
00016 const Teuchos::RCP<Albany::Layouts>& dl) :
00017 coordVec("Coord Vec", dl->qp_vector),
00018 weights("Weights", dl->qp_scalar)
00019 {
00020
00021 Teuchos::ParameterList* plist =
00022 p.get<Teuchos::ParameterList*>("Parameter List");
00023 Teuchos::RCP<const Teuchos::ParameterList> reflist =
00024 this->getValidResponseParameters();
00025 plist->validateParameters(*reflist,0);
00026
00027
00028 Teuchos::RCP<PHX::DataLayout> scalar_dl = dl->qp_scalar;
00029 Teuchos::RCP<PHX::DataLayout> vector_dl = dl->qp_vector;
00030
00031 std::vector<PHX::DataLayout::size_type> dims;
00032 vector_dl->dimensions(dims);
00033 numQPs = dims[1];
00034 numDims = dims[2];
00035
00037 Teuchos::RCP<QCAD::MaterialDatabase> materialDB;
00038 Teuchos::RCP<Teuchos::ParameterList> paramsFromProblem =
00039 p.get< Teuchos::RCP<Teuchos::ParameterList> >("Parameters From Problem");
00040 if(paramsFromProblem != Teuchos::null)
00041 materialDB = paramsFromProblem->get< Teuchos::RCP<QCAD::MaterialDatabase> >("MaterialDB");
00042 else materialDB = Teuchos::null;
00043
00044
00045 fieldName = plist->get<std::string>("Field Name");
00046 opRegion = Teuchos::rcp( new QCAD::MeshRegion<EvalT, Traits>("Coord Vec","Weights",*plist,materialDB,dl) );
00047
00048
00049 PHX::MDField<ScalarT> f(fieldName, scalar_dl); field = f;
00050
00051
00052 this->addDependentField(field);
00053 this->addDependentField(coordVec);
00054 this->addDependentField(weights);
00055 opRegion->addDependentFields(this);
00056
00057 this->setName(fieldName+" Response Field Average"+PHX::TypeString<EvalT>::value);
00058
00059 using PHX::MDALayout;
00060
00061
00062 p.set("Stand-alone Evaluator", false);
00063 std::string local_response_name =
00064 fieldName + " Local Response Field Average";
00065 std::string global_response_name =
00066 fieldName + " Global Response Field Average";
00067 int worksetSize = scalar_dl->dimension(0);
00068 int responseSize = 2;
00069 Teuchos::RCP<PHX::DataLayout> local_response_layout =
00070 Teuchos::rcp(new MDALayout<Cell,Dim>(worksetSize, responseSize));
00071 Teuchos::RCP<PHX::DataLayout> global_response_layout =
00072 Teuchos::rcp(new MDALayout<Dim>(responseSize));
00073 PHX::Tag<ScalarT> local_response_tag(local_response_name,
00074 local_response_layout);
00075 PHX::Tag<ScalarT> global_response_tag(global_response_name,
00076 global_response_layout);
00077 p.set("Local Response Field Tag", local_response_tag);
00078 p.set("Global Response Field Tag", global_response_tag);
00079 PHAL::SeparableScatterScalarResponse<EvalT,Traits>::setup(p,dl);
00080 }
00081
00082
00083 template<typename EvalT, typename Traits>
00084 void QCAD::ResponseFieldAverage<EvalT, Traits>::
00085 postRegistrationSetup(typename Traits::SetupData d,
00086 PHX::FieldManager<Traits>& fm)
00087 {
00088 this->utils.setFieldData(field,fm);
00089 this->utils.setFieldData(coordVec,fm);
00090 this->utils.setFieldData(weights,fm);
00091 opRegion->postRegistrationSetup(fm);
00092 PHAL::SeparableScatterScalarResponse<EvalT,Traits>::postRegistrationSetup(d,fm);
00093 }
00094
00095
00096 template<typename EvalT, typename Traits>
00097 void QCAD::ResponseFieldAverage<EvalT, Traits>::
00098 preEvaluate(typename Traits::PreEvalData workset)
00099 {
00100 for (typename PHX::MDField<ScalarT>::size_type i=0;
00101 i<this->global_response.size(); i++)
00102 this->global_response[i] = 0.0;
00103
00104
00105 PHAL::SeparableScatterScalarResponse<EvalT,Traits>::preEvaluate(workset);
00106 }
00107
00108
00109 template<typename EvalT, typename Traits>
00110 void QCAD::ResponseFieldAverage<EvalT, Traits>::
00111 evaluateFields(typename Traits::EvalData workset)
00112 {
00113
00114 for (typename PHX::MDField<ScalarT>::size_type i=0;
00115 i<this->local_response.size(); i++)
00116 this->local_response[i] = 0.0;
00117
00118 ScalarT t;
00119
00120 if(!opRegion->elementBlockIsInRegion(workset.EBName))
00121 return;
00122
00123 for (std::size_t cell=0; cell < workset.numCells; ++cell)
00124 {
00125 if(!opRegion->cellIsInRegion(cell)) continue;
00126
00127
00128 for (std::size_t qp=0; qp < numQPs; ++qp) {
00129 t = field(cell,qp) * weights(cell,qp);
00130 this->local_response(cell,0) += t;
00131 this->global_response(0) += t;
00132
00133 t = weights(cell,qp);
00134 this->local_response(cell,1) += t;
00135 this->global_response(1) += t;
00136 }
00137
00138 }
00139
00140
00141 PHAL::SeparableScatterScalarResponse<EvalT,Traits>::evaluateFields(workset);
00142 }
00143
00144
00145 template<typename EvalT, typename Traits>
00146 void QCAD::ResponseFieldAverage<EvalT, Traits>::
00147 postEvaluate(typename Traits::PostEvalData workset)
00148 {
00149
00150 Teuchos::RCP< Teuchos::ValueTypeSerializer<int,ScalarT> > serializer =
00151 workset.serializerManager.template getValue<EvalT>();
00152 Teuchos::reduceAll(
00153 *workset.comm, *serializer, Teuchos::REDUCE_SUM,
00154 this->global_response.size(), &this->global_response[0],
00155 &this->global_response[0]);
00156
00157 int iNormalizer = 1;
00158 if( fabs(this->global_response[iNormalizer]) > 1e-9 ) {
00159 for( int i=0; i < this->global_response.size(); i++) {
00160 if( i == iNormalizer ) continue;
00161 this->global_response[i] /= this->global_response[iNormalizer];
00162 }
00163
00164
00165
00166 }
00167
00168
00169 PHAL::SeparableScatterScalarResponse<EvalT,Traits>::postEvaluate(workset);
00170 }
00171
00172
00173 template<typename EvalT,typename Traits>
00174 Teuchos::RCP<const Teuchos::ParameterList>
00175 QCAD::ResponseFieldAverage<EvalT,Traits>::getValidResponseParameters() const
00176 {
00177 Teuchos::RCP<Teuchos::ParameterList> validPL =
00178 rcp(new Teuchos::ParameterList("Valid ResponseFieldAverage Params"));;
00179 Teuchos::RCP<const Teuchos::ParameterList> baseValidPL =
00180 PHAL::SeparableScatterScalarResponse<EvalT,Traits>::getValidResponseParameters();
00181 validPL->setParameters(*baseValidPL);
00182
00183 Teuchos::RCP<const Teuchos::ParameterList> regionValidPL =
00184 QCAD::MeshRegion<EvalT,Traits>::getValidParameters();
00185 validPL->setParameters(*regionValidPL);
00186
00187 validPL->set<std::string>("Name", "", "Name of response function");
00188 validPL->set<int>("Phalanx Graph Visualization Detail", 0, "Make dot file to visualize phalanx graph");
00189 validPL->set<std::string>("Field Name", "", "Scalar field to average");
00190 validPL->set<std::string>("Description", "", "Description of this response used by post processors");
00191
00192 return validPL;
00193 }
00194
00195
00196