Split BA81 into a fit function and expectation
[openmx:openmx.git] / src / omxMLFitFunction.cpp
1  /*
2  *  Copyright 2007-2013 The OpenMx Project
3  *
4  *  Licensed under the Apache License, Version 2.0 (the "License");
5  *  you may not use this file except in compliance with the License.
6  *  You may obtain a copy of the License at
7  *
8  *       http://www.apache.org/licenses/LICENSE-2.0
9  *
10  *   Unless required by applicable law or agreed to in writing, software
11  *   distributed under the License is distributed on an "AS IS" BASIS,
12  *   WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13  *  See the License for the specific language governing permissions and
14  *  limitations under the License.
15  */
16
17 #include <R.h>
18 #include <Rinternals.h>
19 #include <Rdefines.h>
20 #include <R_ext/Rdynload.h>
21 #include <R_ext/BLAS.h>
22 #include <R_ext/Lapack.h>
23 #include "omxAlgebraFunctions.h"
24 #include "omxExpectation.h"
25 #include "omxMLFitFunction.h"
26 #include "omxFIMLFitFunction.h"
27 #include "omxRAMExpectation.h"
28
29 void omxDestroyMLFitFunction(omxFitFunction *oo) {
30
31         if(OMX_DEBUG) {mxLog("Freeing ML Fit Function.");}
32         omxMLFitFunction* omlo = ((omxMLFitFunction*)oo->argStruct);
33
34         if(omlo->localCov != NULL)      omxFreeMatrixData(omlo->localCov);
35         if(omlo->localProd != NULL)     omxFreeMatrixData(omlo->localProd);
36         if(omlo->P != NULL)                     omxFreeMatrixData(omlo->P);
37         if(omlo->C != NULL)                     omxFreeMatrixData(omlo->C);
38         if(omlo->I != NULL)                     omxFreeMatrixData(omlo->I);
39 }
40
41 omxRListElement* omxSetFinalReturnsMLFitFunction(omxFitFunction *oo, int *numReturns) {
42         *numReturns = 3;
43         omxRListElement* retVal = (omxRListElement*) R_alloc(*numReturns, sizeof(omxRListElement));
44         double det = 0.0;
45         omxMatrix* cov = ((omxMLFitFunction*)oo->argStruct)->observedCov;
46         int ncols = ((omxMLFitFunction*)oo->argStruct)->observedCov->cols;
47     
48         retVal[0].numValues = 1;
49         retVal[0].values = (double*) R_alloc(1, sizeof(double));
50         strncpy(retVal[0].label, "Minus2LogLikelihood", 20);
51         retVal[0].values[0] = omxMatrixElement(oo->matrix, 0, 0);
52
53         retVal[1].numValues = 1;
54         retVal[1].values = (double*) R_alloc(1, sizeof(double));
55         strncpy(retVal[1].label, "SaturatedLikelihood", 20);
56         retVal[1].values[0] = (((omxMLFitFunction*)oo->argStruct)->logDetObserved + ncols) * (((omxMLFitFunction*)oo->argStruct)->n - 1);
57
58         retVal[2].numValues = 1;
59         retVal[2].values = (double*) R_alloc(1, sizeof(double));
60         strncpy(retVal[2].label, "IndependenceLikelihood", 23);
61         // Independence model assumes all-zero manifest covariances.
62         // (det(expected) + tr(observed * expected^-1)) * (n - 1);
63         // expected is the diagonal of the observed.  Inverse expected is 1/each diagonal value.
64         // Therefore the diagonal elements of observed * expected^-1 are each 1.
65         // So the trace of the matrix is the same as the number of columns.
66         // The determinant of a diagonal matrix is the product of the diagonal elements.
67         // Since these are the same in the expected as in the observed, we can get 'em direct.
68         for(int i = 0; i < ncols; i++) {
69                 // We sum logs instead of logging the product.
70                 det += log(omxMatrixElement(cov, i, i));
71         }
72         if(OMX_DEBUG) { mxLog("det: %f, tr: %f, n= %d, total:%f", det, ncols, ((omxMLFitFunction*)oo->argStruct)->n, (ncols + det) * (((omxMLFitFunction*)oo->argStruct)->n - 1)); }
73         if(OMX_DEBUG) { omxPrint(cov, "Observed:"); }
74         retVal[2].values[0] = (ncols + det) * (((omxMLFitFunction*)oo->argStruct)->n - 1);
75
76         return retVal;
77 }
78
79 static void omxCallMLFitFunction(omxFitFunction *oo, int want, FitContext *) {
80
81         if (want & FF_COMPUTE_PREOPTIMIZE) return;
82
83         if(OMX_DEBUG) { mxLog("Beginning ML Evaluation.");}
84         // Requires: Data, means, covariances.
85
86         double sum = 0.0, det = 0.0;
87         char u = 'U';
88         char r = 'R';
89         int info = 0;
90         double oned = 1.0;
91         double zerod = 0.0;
92         double minusoned = -1.0;
93         int onei = 1;
94         double fmean = 0.0;
95
96         omxMatrix *scov, *smeans, *cov, *means, *localCov, *localProd, *P, *C;
97
98         omxMLFitFunction *omo = ((omxMLFitFunction*)oo->argStruct);
99         
100     /* Locals for readability.  Compiler should cut through this. */
101         scov            = omo->observedCov;
102         smeans          = omo->observedMeans;
103         cov                     = omo->expectedCov;
104         means           = omo->expectedMeans;
105         localCov        = omo->localCov;
106         localProd       = omo->localProd;
107         P                       = omo->P;
108         C                       = omo->C;
109         double n        = omo->n;
110         double Q        = omo->logDetObserved;
111         omxExpectation* expectation = oo->expectation;
112
113     /* Recompute and recopy */
114         omxExpectationCompute(expectation, NULL);
115
116         omxCopyMatrix(localCov, cov);                           // But expected cov is destroyed in inversion
117
118         if(OMX_DEBUG_ALGEBRA) {
119                 omxPrint(scov, "Observed Covariance is");
120                 omxPrint(localCov, "Implied Covariance Is");
121                 omxPrint(cov, "Original Covariance Is");
122         }
123
124         /* Calculate |expected| */
125
126 //      F77_CALL(dgetrf)(&(localCov->cols), &(localCov->rows), localCov->data, &(localCov->cols), ipiv, &info);
127         F77_CALL(dpotrf)(&u, &(localCov->cols), localCov->data, &(localCov->cols), &info);
128
129         if(OMX_DEBUG_ALGEBRA) { mxLog("Info on LU Decomp: %d", info);}
130         if(info > 0) {
131                 omxRaiseErrorf(oo->matrix->currentState,
132                                "Expected covariance matrix is non-positive-definite after %d evaluations",
133                                oo->matrix->currentState->computeCount);
134                 return;
135         }
136
137         //det = log(det)        // TVO: changed multiplication of det to first log and the summing up; this line should just set det to zero.
138         for(info = 0; info < localCov->cols; info++) {          // |cov| is the square of the product of the diagonal elements of U from the LU factorization.
139                 det += log(fabs(localCov->data[info+localCov->rows*info])); // TVO: changed * to + and added fabs command
140         }
141         det *= 2.0;             // TVO: instead of det *= det;
142
143         if(OMX_DEBUG_ALGEBRA) { mxLog("Determinant of Expected Cov: %f", exp(det)); }
144         // TVO: removed det = log(fabs(det))
145         if(OMX_DEBUG_ALGEBRA) { mxLog("Log of Determinant of Expected Cov: %f", det); }
146
147         /* Calculate Expected^(-1) */
148 //      F77_CALL(dgetri)(&(localCov->rows), localCov->data, &(localCov->cols), ipiv, work, lwork, &info);
149         F77_CALL(dpotri)(&u, &(localCov->rows), localCov->data, &(localCov->cols), &info);
150         if(OMX_DEBUG_ALGEBRA) { mxLog("Info on Invert: %d", info); }
151
152         if(OMX_DEBUG_ALGEBRA) {omxPrint(cov, "Expected Covariance Matrix:");}
153         if(OMX_DEBUG_ALGEBRA) {omxPrint(localCov, "Inverted Matrix:");}
154
155         /* Calculate C = Observed * expected^(-1) */
156
157         if(OMX_DEBUG_ALGEBRA) {mxLog("Call is: DSYMM(%d, %d, %f, %0x, %d, %0x, %d, %f, %0x, %d)",
158                                         (scov->rows), (localCov->cols), oned, scov->data, (localCov->leading),
159                                         localCov->data, (localCov->leading), zerod, localProd->data, (localProd->leading));}
160
161
162         // Stop gcc from issuing a warning
163         int majority = *(scov->majority) == 'n' ? scov->rows : scov->cols;
164
165         /*  TODO:  Make sure leading edges are being appropriately calculated, and sub them back into this */
166         F77_CALL(dsymm)(&r, &u, &(localCov->rows), &(scov->cols),
167                                         &oned, localCov->data, &(majority),
168                                         scov->data, &(majority),
169                                         &zerod, localProd->data, &(localProd->leading));
170
171     /* And get the trace of the result */
172
173         for(info = 0; info < localCov->cols; info++) {
174                 sum += localProd->data[info*localCov->cols + info];
175         }
176
177 //      for(info = 0; info < (localCov->cols * localCov->rows); info++) {
178 //              sum += localCov->data[info] * scov->data[info];
179 //      }
180
181         if(OMX_DEBUG_ALGEBRA) {omxPrint(scov, "Observed Covariance Matrix:");}
182         if(OMX_DEBUG_ALGEBRA) {omxPrint(localCov, "Inverse Matrix:");}
183         if(OMX_DEBUG_ALGEBRA) {omxPrint(localProd, "Product Matrix:");}
184
185         if(means != NULL) {
186                 if(OMX_DEBUG_ALGEBRA) { mxLog("Means Likelihood Calculation"); }
187                 omxRecompute(means);
188                 omxCopyMatrix(P, means);
189                 // P = means - smeans
190                 if(OMX_DEBUG_ALGEBRA) {omxPrint(means, "means");}
191                 if(OMX_DEBUG_ALGEBRA) {omxPrint(smeans, "smeans");}
192                 F77_CALL(daxpy)(&(smeans->cols), &minusoned, smeans->data, &onei, P->data, &onei);
193                 if(OMX_DEBUG_ALGEBRA) {omxPrint(P, "means - smeans");}
194                 // C = P * Cov
195                 F77_CALL(dsymv)(&u, &(localCov->rows), &oned, localCov->data, &(localCov->leading), P->data, &onei, &zerod, C->data, &onei);
196                 // P = C * P'
197                 fmean = F77_CALL(ddot)(&(C->cols), P->data, &onei, C->data, &onei);
198
199                 if(OMX_DEBUG_ALGEBRA) { mxLog("Mean contribution to likelihood is %f per row.", fmean); }
200                 if(fmean < 0.0) fmean = 0.0;
201         }
202
203         oo->matrix->data[0] = (sum + det) * (n - 1) + fmean * (n);
204
205         if(OMX_DEBUG) { mxLog("MLFitFunction value comes to: %f (Chisq: %f).", oo->matrix->data[0], (sum + det) - Q - cov->cols); }
206
207 }
208
209 void omxPopulateMLAttributes(omxFitFunction *oo, SEXP algebra) {
210     if(OMX_DEBUG) { mxLog("Populating ML Attributes."); }
211
212         omxMLFitFunction *argStruct = ((omxMLFitFunction*)oo->argStruct);
213         omxMatrix *expCovInt = argStruct->expectedCov;                  // Expected covariance
214         omxMatrix *expMeanInt = argStruct->expectedMeans;                       // Expected means
215
216         SEXP expCovExt, expMeanExt, gradients;
217         PROTECT(expCovExt = allocMatrix(REALSXP, expCovInt->rows, expCovInt->cols));
218         for(int row = 0; row < expCovInt->rows; row++)
219                 for(int col = 0; col < expCovInt->cols; col++)
220                         REAL(expCovExt)[col * expCovInt->rows + row] =
221                                 omxMatrixElement(expCovInt, row, col);
222
223         if (expMeanInt != NULL) {
224                 PROTECT(expMeanExt = allocMatrix(REALSXP, expMeanInt->rows, expMeanInt->cols));
225                 for(int row = 0; row < expMeanInt->rows; row++)
226                         for(int col = 0; col < expMeanInt->cols; col++)
227                                 REAL(expMeanExt)[col * expMeanInt->rows + row] =
228                                         omxMatrixElement(expMeanInt, row, col);
229         } else {
230                 PROTECT(expMeanExt = allocMatrix(REALSXP, 0, 0));               
231         }   
232
233         if (0) {
234                 // TODO, fix for new gradient internal API
235         } else {
236                 PROTECT(gradients = allocMatrix(REALSXP, 0, 0));
237         }
238     
239         setAttrib(algebra, install("expCov"), expCovExt);
240         setAttrib(algebra, install("expMean"), expMeanExt);
241         setAttrib(algebra, install("gradients"), gradients);
242         
243         UNPROTECT(3);
244
245 }
246
247 void omxSetMLFitFunctionCalls(omxFitFunction* oo) {
248         
249         /* Set FitFunction Calls to ML FitFunction Calls */
250         oo->computeFun = omxCallMLFitFunction;
251         oo->destructFun = omxDestroyMLFitFunction;
252         oo->setFinalReturns = omxSetFinalReturnsMLFitFunction;
253         oo->populateAttrFun = omxPopulateMLAttributes;
254 }
255
256
257 void omxInitMLFitFunction(omxFitFunction* oo) {
258
259         if(OMX_DEBUG) { mxLog("Initializing ML fit function."); }
260
261         int info = 0;
262         double det = 1.0;
263         char u = 'U';
264         
265         /* Read and set expectation */
266         omxSetMLFitFunctionCalls(oo);
267
268         if (!oo->expectation) { error("%s requires an expectation", oo->fitType); }
269
270         omxData* dataMat = oo->expectation->data;
271
272         if(!(dataMat == NULL) && strncmp(omxDataType(dataMat), "cov", 3) != 0 && strncmp(omxDataType(dataMat), "cor", 3) != 0) {
273                 if(strncmp(omxDataType(dataMat), "raw", 3) == 0) {
274                         if(OMX_DEBUG) { mxLog("Raw Data: Converting to FIML."); }
275                         omxInitFIMLFitFunction(oo);
276                         return;
277                 }
278                 char *errstr = (char*) calloc(250, sizeof(char));
279                 sprintf(errstr, "ML FitFunction unable to handle data type %s.\n", omxDataType(dataMat));
280                 omxRaiseError(oo->matrix->currentState, -1, errstr);
281                 free(errstr);
282                 if(OMX_DEBUG) { mxLog("ML FitFunction unable to handle data type %s.  Aborting.", omxDataType(dataMat)); }
283                 return;
284         }
285
286         omxMLFitFunction *newObj = (omxMLFitFunction*) R_alloc(1, sizeof(omxMLFitFunction));
287         oo->argStruct = (void*)newObj;
288
289         if(OMX_DEBUG) { mxLog("Processing Observed Covariance."); }
290         newObj->observedCov = omxDataMatrix(dataMat, NULL);
291         if(OMX_DEBUG) { mxLog("Processing Observed Means."); }
292         newObj->observedMeans = omxDataMeans(dataMat, NULL, NULL);
293         if(OMX_DEBUG && newObj->observedMeans == NULL) { mxLog("ML: No Observed Means."); }
294         if(OMX_DEBUG) { mxLog("Processing n."); }
295         newObj->n = omxDataNumObs(dataMat);
296
297         newObj->expectedCov = omxGetExpectationComponent(oo->expectation, oo, "cov");
298         newObj->expectedMeans = omxGetExpectationComponent(oo->expectation, oo, "means");
299
300         if(newObj->expectedCov == NULL) {
301                 omxRaiseError(oo->matrix->currentState, OMX_DEVELOPER_ERROR,
302                         "Developer Error in ML-based fit function object: ML's expectation must specify a model-implied covariance matrix.\nIf you are not developing a new expectation type, you should probably post this to the OpenMx forums.");
303                 return;
304         }
305
306         // Error Checking: Observed/Expected means must agree.  
307         // ^ is XOR: true when one is false and the other is not.
308         if((newObj->expectedMeans == NULL) ^ (newObj->observedMeans == NULL)) {
309                 if(newObj->expectedMeans != NULL) {
310                         omxRaiseError(oo->matrix->currentState, OMX_ERROR,
311                                 "Observed means not detected, but an expected means matrix was specified.\n  If you provide observed means, you must specify a model for the means.\n");
312                         return;
313                 } else {
314                         omxRaiseError(oo->matrix->currentState, OMX_ERROR,
315                                 "Observed means were provided, but an expected means matrix was not specified.\n  If you  wish to model the means, you must provide observed means.\n");
316                         return;         
317                 }
318         }
319
320         /* Temporary storage for calculation */
321         int rows = newObj->observedCov->rows;
322         int cols = newObj->observedCov->cols;
323         newObj->localCov = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
324         newObj->localProd = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
325         newObj->P = omxInitMatrix(NULL, 1, cols, TRUE, oo->matrix->currentState);
326         newObj->C = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
327         newObj->I = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
328
329         for(int i = 0; i < rows; i++) omxSetMatrixElement(newObj->I, i, i, 1.0);
330
331         omxCopyMatrix(newObj->localCov, newObj->observedCov);
332
333         newObj->lwork = newObj->expectedCov->rows;
334         newObj->work = (double*)R_alloc(newObj->lwork, sizeof(double));
335
336         // TODO: Determine where the saturated model computation should go.
337
338         F77_CALL(dpotrf)(&u, &(newObj->localCov->cols), newObj->localCov->data, &(newObj->localCov->cols), &info);
339
340         if(OMX_DEBUG) { mxLog("Info on LU Decomp: %d", info); }
341         if(info != 0) {
342                 char *errstr = (char*) calloc(250, sizeof(char));
343                 sprintf(errstr, "Observed Covariance Matrix is non-positive-definite.\n");
344                 omxRaiseError(oo->matrix->currentState, -1, errstr);
345                 free(errstr);
346                 return;
347         }
348         for(info = 0; info < newObj->localCov->cols; info++) {
349                 det *= omxMatrixElement(newObj->localCov, info, info);
350         }
351         det *= det;                                     // Product of squares.
352
353         if(OMX_DEBUG) { mxLog("Determinant of Observed Cov: %f", det); }
354         newObj->logDetObserved = log(det);
355         if(OMX_DEBUG) { mxLog("Log Determinant of Observed Cov: %f", newObj->logDetObserved); }
356
357         omxCopyMatrix(newObj->localCov, newObj->expectedCov);
358 }
359
360 void omxSetMLFitFunctionGradient(omxFitFunction* oo, void (*derivativeFun)(omxFitFunction*, double*)) {
361     if(strncmp("omxMLFitFunction", oo->fitType, 16)) {
362         char errorstr[250];
363         sprintf(errorstr, "PROGRAMMER ERROR: Using vanilla-ML gradient with FitFunction of type %s", oo->fitType);
364         omxRaiseError(oo->matrix->currentState, -2, errorstr);
365         return;
366     }
367     
368     if(derivativeFun == NULL) {
369         char errorstr[250];
370         sprintf(errorstr, "Programmer error: ML gradient given NULL gradient function.");
371         omxRaiseError(oo->matrix->currentState, -2, errorstr);
372         return;
373     }
374     
375     //oo->gradientFun = derivativeFun; TODO
376 }
377
378 void omxSetMLFitFunctionGradientComponents(omxFitFunction* oo, void (*derivativeFun)(omxFitFunction*, omxMatrix**, omxMatrix**, int*)) {
379     if(OMX_DEBUG) { mxLog("Setting up gradient component function for ML FitFunction."); }
380     if(!strncmp("omxFIMLFitFunction", oo->fitType, 16)) {
381         if(OMX_DEBUG) { mxLog("FIML FitFunction gradients not yet implemented. Skipping."); }
382         return; // ERROR:NYI.
383     }
384     
385     if(derivativeFun == NULL) {
386         char errorstr[250];
387         sprintf(errorstr, "Programmer error: ML gradient components given NULL gradient function.");
388         omxRaiseError(oo->matrix->currentState, -2, errorstr);
389         return;
390     }
391     
392     omxMLFitFunction *omo = ((omxMLFitFunction*) oo->argStruct);
393     int rows = omo->observedCov->rows;
394     int cols = omo->observedCov->cols;
395     size_t nFreeVars = oo->freeVarGroup->vars.size();
396             
397     omo->derivativeFun = derivativeFun;
398     omo->X  = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
399     omo->Y  = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
400     omo->Ms = omxInitMatrix(NULL, 1, cols, TRUE, oo->matrix->currentState);
401     omo->Mu = omxInitMatrix(NULL, 1, cols, TRUE, oo->matrix->currentState);
402     omo->dSigma = (omxMatrix**) R_alloc(nFreeVars, sizeof(omxMatrix*));
403     omo->dMu = (omxMatrix**) R_alloc(nFreeVars, sizeof(omxMatrix*));
404     for(size_t i = 0; i < nFreeVars; i++) {
405         omo->dSigma[i] = omxInitMatrix(NULL, rows, cols, TRUE, oo->matrix->currentState);
406         omo->dMu[i] = omxInitMatrix(NULL, rows, 1, TRUE, oo->matrix->currentState);
407     }
408     //oo->gradientFun = omxCalculateMLGradient; TODO
409 }
410
411 void omxCalculateMLGradient(omxFitFunction* oo, double* gradient) {
412
413     if(OMX_DEBUG) { mxLog("Beginning ML Gradient Calculation."); }
414     // mxLog("Beginning ML Gradient Calculation, Iteration %d.%d (%d)\n", 
415         // oo->matrix->currentState->majorIteration, oo->matrix->currentState->minorIteration,
416         // oo->matrix->currentState->computeCount); //:::DEBUG:::
417     // 1) Calculate current Expected Covariance
418     // 2) Calculate eCov, the Inverse Expected Covariance matrix 
419     // 3) Calculate C = I - eCov D, where D is the observed covariance matrix
420     // 4) Calculate b = M - [observed M]
421     // 5) For each location in locs:
422     //   gradient[loc] = tr(eCov^-1 %*% dEdt %*% C) - (b^T %*% eCov^-1 %*% dEdt + 2 dMdt^T))eCov^-1 b)
423
424     omxMLFitFunction *omo = ((omxMLFitFunction*)oo->argStruct);
425     
426     /* Locals for readability.  Compiler should cut through this. */
427     omxMatrix *scov         = omo->observedCov;
428     omxMatrix *smeans       = omo->observedMeans;
429     omxMatrix *cov          = omo->expectedCov;
430     omxMatrix *M            = omo->expectedMeans;
431     omxMatrix *eCov         = omo->localCov;        // TODO: Maybe need to avoid reusing these
432     omxMatrix *I            = omo->I;
433     omxMatrix *C            = omo->C;
434     omxMatrix *X            = omo->X;
435     omxMatrix *Y            = omo->Y;
436     omxMatrix *Mu           = omo->Mu;
437     omxMatrix *Ms           = omo->Ms;
438     omxMatrix *P            = omo->P;
439     double n                = omo->n;
440     omxMatrix** dSigmas     = omo->dSigma;
441     omxMatrix** dMus        = omo->dMu;
442     
443     size_t gradientSize = oo->freeVarGroup->vars.size();
444     
445     char u = 'U';
446     int info;
447     double minusoned = -1.0;
448     int onei = 1;
449     int status[gradientSize];
450     int nLocs = gradientSize;
451     
452     // Calculate current FitFunction values
453     // We can safely assume this has been done
454     // omxFitFunctionCompute(oo);
455     
456     // Calculate current eCov
457     
458     omxCopyMatrix(eCov, cov);                           // But expected cov is destroyed in inversion
459     
460     F77_CALL(dpotrf)(&u, &(eCov->cols), eCov->data, &(eCov->cols), &info);
461
462     if(OMX_DEBUG_ALGEBRA) { mxLog("Info on LU Decomp: %d", info);}
463     if(info > 0) {
464             char *errstr = (char*) calloc(250, sizeof(char));
465         sprintf(errstr, "Expected covariance matrix is non-positive-definite");
466         if(oo->matrix->currentState->computeCount <= 0) {
467             strncat(errstr, " at starting values", 20);
468         }
469         strncat(errstr, ".\n", 3);
470         omxRaiseError(oo->matrix->currentState, -1, errstr);                        // Raise error
471         free(errstr);
472         return;                                                                     // Leave output untouched
473     }
474     
475     F77_CALL(dpotri)(&u, &(eCov->rows), eCov->data, &(eCov->cols), &info);
476     if(info > 0) {
477             char *errstr = (char*) calloc(250, sizeof(char));
478         sprintf(errstr, "Expected covariance matrix is not invertible");
479         if(oo->matrix->currentState->computeCount <= 0) {
480             strncat(errstr, " at starting values", 20);
481         }
482         strncat(errstr, ".\n", 3);
483         omxRaiseError(oo->matrix->currentState, -1, errstr);                        // Raise error
484         free(errstr);
485         return;
486     }
487     // Calculate P = expected means - observed means
488     if(M != NULL) {
489         omxCopyMatrix(P, smeans);
490         F77_CALL(daxpy)(&(smeans->cols), &minusoned, M->data, &onei, P->data, &onei);
491     }
492         
493         // Reset C and Calculate C = I - eCov * oCov
494     omxCopyMatrix(C, I);
495     omxDSYMM(TRUE, -1.0, eCov, scov, 1.0, C);
496     
497     // For means, calculate Ms = eCov-1 %*% P
498     if(M != NULL)
499         omxDSYMM(FALSE, 1.0, eCov, P, 0.0, Ms);
500     
501     // Calculate parameter-level derivatives
502     // TODO: Parallelize Here.
503
504     if(OMX_DEBUG)  { mxLog("Calling component function."); }
505     omo->derivativeFun(oo, dSigmas, dMus, status);
506     
507     for(int currentLoc = 0; currentLoc < nLocs; currentLoc++) {
508         double meanInfluence, covInfluence;
509         if(status[currentLoc] < 0) continue;  // Failure in computation--skip.
510         //   gradient[loc] = tr(eCov^-1 %*% dEdt %*% C) - 
511         //    (b^T %*% eCov^-1 %*% dEdt + 2 dMdt^T))eCov^-1 b)
512         // omxDGEMM(FALSE, FALSE, 1.0, dSigmas[currentLoc], C, 0.0, Y);
513         omxDSYMM(TRUE, 1.0, dSigmas[currentLoc], C, 0.0, Y);
514         omxDSYMM(TRUE, 1.0, eCov, Y, 0.0, X);
515         gradient[currentLoc] = 0;
516         covInfluence = 0.0;
517         for(int i = 0; i < eCov->cols; i++) 
518             covInfluence += omxMatrixElement(X, i, i);
519         if(M != NULL) {
520             omxCopyMatrix(Mu, dMus[currentLoc]);
521             omxDSYMV(1.0, dSigmas[currentLoc], Ms, 2.0, Mu);
522             meanInfluence = F77_CALL(ddot)(&(eCov->cols), Mu->data, &onei, Ms->data, &onei);
523         } else {
524             meanInfluence = 0;
525         }
526         gradient[currentLoc] = (covInfluence * (n-1)) - (meanInfluence * n);
527         if(OMX_DEBUG) { 
528             mxLog("Calculation for Gradient value %d: Cov: %3.9f, Mean: %3.9f, total: %3.9f",
529             currentLoc, covInfluence, meanInfluence, gradient[currentLoc]); 
530         }
531     }
532 }